File indexing completed on 2024-04-06 12:31:58
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0012
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017
0018 #include "DataFormats/CTPPSDetId/interface/CTPPSDetId.h"
0019 #include "DataFormats/CTPPSReco/interface/CTPPSLocalTrackLite.h"
0020 #include "DataFormats/CTPPSReco/interface/CTPPSLocalTrackLiteFwd.h"
0021
0022 #include "TFile.h"
0023 #include "TH1D.h"
0024 #include "TH2D.h"
0025
0026 #include <map>
0027
0028
0029
0030 class CTPPSDirectProtonSimulationValidator : public edm::one::EDAnalyzer<> {
0031 public:
0032 explicit CTPPSDirectProtonSimulationValidator(const edm::ParameterSet&);
0033
0034 private:
0035 void beginJob() override {}
0036 void analyze(const edm::Event&, const edm::EventSetup&) override;
0037 void endJob() override;
0038
0039 edm::EDGetTokenT<CTPPSLocalTrackLiteCollection> simuTracksToken_;
0040 edm::EDGetTokenT<CTPPSLocalTrackLiteCollection> recoTracksToken_;
0041
0042 std::string outputFile_;
0043
0044 struct RPPlots {
0045 std::unique_ptr<TH2D> h2_xr_vs_xs, h2_yr_vs_ys;
0046 std::unique_ptr<TH1D> h_de_x, h_de_y;
0047
0048 RPPlots()
0049 : h2_xr_vs_xs(new TH2D("", ";x_simu (mm);x_reco (mm)", 100, -10., +10., 100, -10, +10.)),
0050 h2_yr_vs_ys(new TH2D("", "y_simu (mm);y_reco (mm)", 100, -10., +10., 100, -10, +10.)),
0051 h_de_x(new TH1D("", ";x (mm)", 200, -100E-3, +100E-3)),
0052 h_de_y(new TH1D("", ";y (mm)", 200, -100E-3, +100E-3)) {}
0053
0054 void fill(double simu_x, double simu_y, double reco_x, double reco_y) {
0055 h2_xr_vs_xs->Fill(simu_x, reco_x);
0056 h2_yr_vs_ys->Fill(simu_y, reco_y);
0057
0058 h_de_x->Fill(reco_x - simu_x);
0059 h_de_y->Fill(reco_y - simu_y);
0060 }
0061
0062 void write() const {
0063 h2_xr_vs_xs->Write("h2_xr_vs_xs");
0064 h2_yr_vs_ys->Write("h2_yr_vs_ys");
0065 h_de_x->Write("h_de_x");
0066 h_de_y->Write("h_de_y");
0067 }
0068 };
0069
0070 std::map<unsigned int, RPPlots> rpPlots_;
0071 };
0072
0073
0074
0075 CTPPSDirectProtonSimulationValidator::CTPPSDirectProtonSimulationValidator(const edm::ParameterSet& iConfig)
0076 : simuTracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("simuTracksTag"))),
0077 recoTracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("recoTracksTag"))),
0078 outputFile_(iConfig.getParameter<std::string>("outputFile")) {}
0079
0080
0081
0082 void CTPPSDirectProtonSimulationValidator::analyze(const edm::Event& iEvent, const edm::EventSetup&) {
0083
0084 edm::Handle<CTPPSLocalTrackLiteCollection> simuTracks;
0085 iEvent.getByToken(simuTracksToken_, simuTracks);
0086
0087 edm::Handle<CTPPSLocalTrackLiteCollection> recoTracks;
0088 iEvent.getByToken(recoTracksToken_, recoTracks);
0089
0090
0091 for (const auto& simuTrack : *simuTracks) {
0092 const CTPPSDetId rpId(simuTrack.rpId());
0093 for (const auto& recoTrack : *recoTracks) {
0094 if (simuTrack.rpId() == recoTrack.rpId()) {
0095 unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
0096 rpPlots_[rpDecId].fill(simuTrack.x(), simuTrack.y(), recoTrack.x(), recoTrack.y());
0097 }
0098 }
0099 }
0100 }
0101
0102
0103
0104 void CTPPSDirectProtonSimulationValidator::endJob() {
0105 auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
0106
0107 for (const auto& it : rpPlots_) {
0108 gDirectory = f_out->mkdir(Form("RP %u", it.first));
0109 it.second.write();
0110 }
0111 }
0112
0113
0114
0115 DEFINE_FWK_MODULE(CTPPSDirectProtonSimulationValidator);