CMS 3D CMS Logo

CTPPSBeamSmearingValidator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
8 
11 
13 
14 #include "TFile.h"
15 #include "TH1D.h"
16 
17 #include <map>
18 
19 //----------------------------------------------------------------------------------------------------
20 
22 public:
24 
26 
27 private:
28  void analyze(const edm::Event &, const edm::EventSetup &) override;
29 
30  void endJob() override;
31 
34 
36 
37  std::unique_ptr<TH1D> h_de_vtx_x_, h_de_vtx_y_, h_de_vtx_z_;
38 
39  struct SectorPlots {
40  std::unique_ptr<TH1D> h_de_th_x, h_de_th_y, h_de_p;
41 
43  : h_de_th_x(new TH1D("", ";#Delta#theta_{x} (rad)", 100, 0., 0.)),
44  h_de_th_y(new TH1D("", ";#Delta#theta_{y} (rad)", 100, 0., 0.)),
45  h_de_p(new TH1D("", ";#Deltap (GeV)", 100, 0., 0.)) {}
46 
47  void write() const {
48  h_de_th_x->Write("h_de_th_x");
49  h_de_th_y->Write("h_de_th_y");
50  h_de_p->Write("h_de_p");
51  }
52  };
53 
54  std::map<unsigned int, SectorPlots> sectorPlots_;
55 };
56 
57 //----------------------------------------------------------------------------------------------------
58 
59 using namespace std;
60 using namespace edm;
61 using namespace HepMC;
62 
63 //----------------------------------------------------------------------------------------------------
64 
66  : tokenBeforeSmearing_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagBeforeSmearing"))),
67  tokenAfterSmearing_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagAfterSmearing"))),
68  outputFile_(iConfig.getParameter<string>("outputFile")),
69  h_de_vtx_x_(new TH1D("h_de_vtx_x", ";#Delta vtx_{x} (mm)", 100, 0., 0.)),
70  h_de_vtx_y_(new TH1D("h_de_vtx_y", ";#Delta vtx_{y} (mm)", 100, 0., 0.)),
71  h_de_vtx_z_(new TH1D("h_de_vtx_z", ";#Delta vtx_{z} (mm)", 100, 0., 0.)) {}
72 
73 //----------------------------------------------------------------------------------------------------
74 
76  // get input
77  edm::Handle<edm::HepMCProduct> hBeforeSmearing;
78  iEvent.getByToken(tokenBeforeSmearing_, hBeforeSmearing);
79  HepMC::GenEvent *orig = (HepMC::GenEvent *)hBeforeSmearing->GetEvent();
80 
81  edm::Handle<edm::HepMCProduct> hAfterSmearing;
82  iEvent.getByToken(tokenAfterSmearing_, hAfterSmearing);
83  HepMC::GenEvent *smear = (HepMC::GenEvent *)hAfterSmearing->GetEvent();
84 
85  // vertices
86  GenEvent::vertex_const_iterator vold, vnew;
87  for (vold = orig->vertices_begin(), vnew = smear->vertices_begin();
88  vold != orig->vertices_end() && vnew != smear->vertices_end();
89  ++vold, ++vnew) {
90  const FourVector &vo = (*vold)->position();
91  const FourVector &vn = (*vnew)->position();
92 
93  // HepMC gives vertex in mm
94  h_de_vtx_x_->Fill(vn.x() - vo.x());
95  h_de_vtx_y_->Fill(vn.y() - vo.y());
96  h_de_vtx_z_->Fill(vn.z() - vo.z());
97  }
98 
99  // particles
100  GenEvent::particle_const_iterator pold, pnew;
101  for (pold = orig->particles_begin(), pnew = smear->particles_begin();
102  pold != orig->particles_end() && pnew != smear->particles_end();
103  ++pold, ++pnew) {
104  FourVector o = (*pold)->momentum(), n = (*pnew)->momentum();
105 
106  // determine direction region
107  signed int idx = -1;
108  const double thetaLim = 0.01; // rad
109  double th = o.theta();
110 
111  if (th < thetaLim)
112  idx = 0;
113  if (th > (M_PI - thetaLim))
114  idx = 1;
115 
116  if (idx < 0)
117  continue;
118 
119  /*
120  cout << "particle\n\told: [" << o.x() << ", " << o.y() << ", " << o.z() << ", " << o.t()
121  << "]\n\tnew: [" << n.x() << ", " << n.y() << ", " << n.z() << ", " << n.t()
122  << "]\n\tregion: " << idx << endl;
123  */
124 
125  // fill histograms
126  auto &sp = sectorPlots_[idx];
127 
128  double othx = o.x() / o.rho(), othy = o.y() / o.rho();
129  double nthx = n.x() / n.rho(), nthy = n.y() / n.rho();
130 
131  sp.h_de_p->Fill(n.rho() - o.rho());
132 
133  sp.h_de_th_x->Fill(nthx - othx);
134  sp.h_de_th_y->Fill(nthy - othy);
135  }
136 }
137 
138 //----------------------------------------------------------------------------------------------------
139 
141  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
142 
143  h_de_vtx_x_->Write();
144  h_de_vtx_y_->Write();
145  h_de_vtx_z_->Write();
146 
147  gDirectory = f_out->mkdir("sector 45");
148  sectorPlots_[0].write();
149 
150  gDirectory = f_out->mkdir("sector 56");
151  sectorPlots_[1].write();
152 }
153 
154 //----------------------------------------------------------------------------------------------------
155 
void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< edm::HepMCProduct > tokenAfterSmearing_
#define M_PI
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
HLT enums.
edm::EDGetTokenT< edm::HepMCProduct > tokenBeforeSmearing_
std::map< unsigned int, SectorPlots > sectorPlots_
CTPPSBeamSmearingValidator(const edm::ParameterSet &)