CMS 3D CMS Logo

CTPPSDirectProtonSimulationValidator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  *
3  * This is a part of CTPPS validation software
4  * Authors:
5  * Jan Kašpar
6  * Laurent Forthomme
7  *
8  ****************************************************************************/
9 
10 
13 
16 
18 
22 
23 #include "TFile.h"
24 #include "TH1D.h"
25 #include "TH2D.h"
26 
27 #include <map>
28 
29 //----------------------------------------------------------------------------------------------------
30 
32 {
33  public:
35 
36  private:
37  void beginJob() override {}
38  void analyze( const edm::Event&, const edm::EventSetup& ) override;
39  void endJob() override;
40 
43 
45 
46  struct RPPlots
47  {
48  std::unique_ptr<TH2D> h2_xr_vs_xs, h2_yr_vs_ys;
49  std::unique_ptr<TH1D> h_de_x, h_de_y;
50 
51  RPPlots() :
52  h2_xr_vs_xs( new TH2D("", "", 100, -10., +10., 100, -10, +10.) ),
53  h2_yr_vs_ys( new TH2D("", "", 100, -10., +10., 100, -10, +10.) ),
54  h_de_x( new TH1D("", "", 100, -0., +0.) ),
55  h_de_y( new TH1D("", "", 100, -0., +0.) )
56  {}
57 
58  void fill(double simu_x, double simu_y, double reco_x, double reco_y)
59  {
60  h2_xr_vs_xs->Fill(simu_x, reco_x);
61  h2_yr_vs_ys->Fill(simu_y, reco_y);
62 
63  h_de_x->Fill(reco_x - simu_x);
64  h_de_y->Fill(reco_y - simu_y);
65  }
66 
67  void write() const
68  {
69  h2_xr_vs_xs->Write("h2_xr_vs_xs");
70  h2_yr_vs_ys->Write("h2_yr_vs_ys");
71  h_de_x->Write("h_de_x");
72  h_de_y->Write("h_de_y");
73  }
74  };
75 
76  std::map<unsigned int, RPPlots> rpPlots_;
77 };
78 
79 //----------------------------------------------------------------------------------------------------
80 
82  simuTracksToken_( consumes<CTPPSLocalTrackLiteCollection>( iConfig.getParameter<edm::InputTag>( "simuTracksTag" ) ) ),
83  recoTracksToken_( consumes<CTPPSLocalTrackLiteCollection>( iConfig.getParameter<edm::InputTag>( "recoTracksTag" ) ) ),
84  outputFile_( iConfig.getParameter<std::string>("outputFile") )
85 {}
86 
87 //----------------------------------------------------------------------------------------------------
88 
90 {
91  // get input
93  iEvent.getByToken( simuTracksToken_, simuTracks );
94 
96  iEvent.getByToken( recoTracksToken_, recoTracks );
97 
98  // process tracks
99  for (const auto& simuTrack : *simuTracks) {
100  const CTPPSDetId rpId(simuTrack.getRPId());
101  for (const auto& recoTrack : *recoTracks) {
102  if (simuTrack.getRPId() == recoTrack.getRPId()) {
103  unsigned int rpDecId = rpId.arm()*100 + rpId.station()*10 + rpId.rp();
104  rpPlots_[rpDecId].fill(simuTrack.getX(), simuTrack.getY(), recoTrack.getX(), recoTrack.getY());
105  }
106  }
107  }
108 }
109 
110 //----------------------------------------------------------------------------------------------------
111 
113 {
114  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
115 
116  for (const auto& it : rpPlots_) {
117  gDirectory = f_out->mkdir(Form("RP %u", it.first));
118  it.second.write();
119  }
120 }
121 
122 //----------------------------------------------------------------------------------------------------
123 
125 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(const edm::Event &, const edm::EventSetup &) override
uint32_t arm() const
Definition: CTPPSDetId.h:51
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > recoTracksToken_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > simuTracksToken_
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
void fill(double simu_x, double simu_y, double reco_x, double reco_y)
HLT enums.