CMS 3D CMS Logo

TrackerHitAnalyzer.cc
Go to the documentation of this file.
2 
9 
12 
14 
15 // tracker info
20 
21 // data in edm::event
24 //#include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
28 
29 // helper files
30 #include "CLHEP/Units/GlobalSystemOfUnits.h"
31 #include <CLHEP/Vector/LorentzVector.h>
32 
33 #include <cstdlib>
34 #include <map>
35 #include <memory>
36 #include <vector>
37 
39  : tGeomEsToken_(esConsumes()),
40  verbose_(ps.getUntrackedParameter<bool>("Verbosity", false)),
41  edmPSimHitContainer_pxlBrlLow_Token_(
42  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("PxlBrlLowSrc"))),
43  edmPSimHitContainer_pxlBrlHigh_Token_(
44  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("PxlBrlHighSrc"))),
45  edmPSimHitContainer_pxlFwdLow_Token_(
46  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("PxlFwdLowSrc"))),
47  edmPSimHitContainer_pxlFwdHigh_Token_(
48  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("PxlFwdHighSrc"))),
49  edmPSimHitContainer_siTIBLow_Token_(
50  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTIBLowSrc"))),
51  edmPSimHitContainer_siTIBHigh_Token_(
52  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTIBHighSrc"))),
53  edmPSimHitContainer_siTOBLow_Token_(
54  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTOBLowSrc"))),
55  edmPSimHitContainer_siTOBHigh_Token_(
56  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTOBHighSrc"))),
57  edmPSimHitContainer_siTIDLow_Token_(
58  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTIDLowSrc"))),
59  edmPSimHitContainer_siTIDHigh_Token_(
60  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTIDHighSrc"))),
61  edmPSimHitContainer_siTECLow_Token_(
62  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTECLowSrc"))),
63  edmPSimHitContainer_siTECHigh_Token_(
64  consumes<edm::PSimHitContainer>(ps.getParameter<edm::InputTag>("SiTECHighSrc"))),
65  edmSimTrackContainerToken_(consumes<edm::SimTrackContainer>(ps.getParameter<edm::InputTag>("G4TrkSrc"))),
66  fDBE(nullptr),
67  conf_(ps),
68  runStandalone(ps.getParameter<bool>("runStandalone")),
69  fOutputFile(ps.getUntrackedParameter<std::string>("outputFile", "TrackerHitHisto.root")),
70  pixelOutput(ps.getParameter<bool>("pixelOutput")) {}
71 
75 
76  Char_t hname1[50], htitle1[80];
77  Char_t hname2[50], htitle2[80];
78  Char_t hname3[50], htitle3[80];
79  Char_t hname4[50], htitle4[80];
80  Char_t hname5[50], htitle5[80];
81  Char_t hname6[50], htitle6[80];
82 
83  if (fDBE != nullptr) {
84  // fDBE->setCurrentFolder("TrackerHitsV/TrackerHitTask");
85 
86  // is there any way to record CPU Info ???
87  // if so, it can be done once - via beginJob()
88  int nbin = 5000;
89 
90  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/");
91  htofeta = ibooker.book2D("tof_eta", "Time of flight vs eta", nbin, -3.0, 3.0, 200, -100, 100);
92  htofphi = ibooker.book2D("tof_phi", "Time of flight vs phi", nbin, -180, 180, 200, -100, 100);
93  htofr = ibooker.book2D("tof_r", "Time of flight vs r", nbin, 0, 300, 200, -100, 100);
94  htofz = ibooker.book2D("tof_z", "Time of flight vs z", nbin, -280, 280, 200, -100, 100);
95 
96  const float E2NEL = 1.;
97 
98  const char *Region[] = {"005", "051", "115", "152", "225", "253", "-050", "-105", "-151", "-215", "-252", "-325"};
99  nbin = 200;
100 
101  // Energy loss histograms
102  for (int i = 0; i < 12; i++) {
103  sprintf(htitle1, "Energy loss in TIB %s", Region[i]);
104  sprintf(htitle2, "Energy loss in TOB %s", Region[i]);
105  sprintf(htitle3, "Energy loss in TID %s", Region[i]);
106  sprintf(htitle4, "Energy loss in TEC %s", Region[i]);
107  if (pixelOutput)
108  sprintf(htitle5, "Energy loss in BPIX %s", Region[i]);
109  if (pixelOutput)
110  sprintf(htitle6, "Energy loss in FPIX %s", Region[i]);
111 
112  sprintf(hname1, "Eloss_TIB_%i", i + 1);
113  sprintf(hname2, "Eloss_TOB_%i", i + 1);
114  sprintf(hname3, "Eloss_TID_%i", i + 1);
115  sprintf(hname4, "Eloss_TEC_%i", i + 1);
116  if (pixelOutput)
117  sprintf(hname5, "Eloss_BPIX_%i", i + 1);
118  if (pixelOutput)
119  sprintf(hname6, "Eloss_FPIX_%i", i + 1);
120 
121  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
122  h1e[i] = ibooker.book1D(hname1, htitle1, nbin, 0.0, 0.001 * E2NEL);
123  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
124  h2e[i] = ibooker.book1D(hname2, htitle2, nbin, 0.0, 0.001 * E2NEL);
125  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
126  h3e[i] = ibooker.book1D(hname3, htitle3, nbin, 0.0, 0.001 * E2NEL);
127  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
128  h4e[i] = ibooker.book1D(hname4, htitle4, nbin, 0.0, 0.001 * E2NEL);
129  if (pixelOutput) {
130  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
131  h5e[i] = ibooker.book1D(hname5, htitle5, nbin, 0.0, 0.001 * E2NEL);
132  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
133  h6e[i] = ibooker.book1D(hname6, htitle6, nbin, 0.0, 0.001 * E2NEL);
134  }
135  }
136 
137  // limits
138  const float high[] = {0.03, 0.03, 0.02, 0.03, 0.03, 0.03};
139  const float low[] = {-0.03, -0.03, -0.02, -0.03, -0.03, -0.03};
140 
141  for (int i = 0; i < 12; i++) {
142  sprintf(htitle1, "Entryx-Exitx in TIB %s", Region[i]);
143  sprintf(htitle2, "Entryx-Exitx in TOB %s", Region[i]);
144  sprintf(htitle3, "Entryx-Exitx in TID %s", Region[i]);
145  sprintf(htitle4, "Entryx-Exitx in TEC %s", Region[i]);
146  if (pixelOutput)
147  sprintf(htitle5, "Entryx-Exitx in BPIX %s", Region[i]);
148  if (pixelOutput)
149  sprintf(htitle6, "Entryx-Exitx in FPIX %s", Region[i]);
150 
151  sprintf(hname1, "Entryx-Exitx_TIB_%i", i + 1);
152  sprintf(hname2, "Entryx-Exitx_TOB_%i", i + 1);
153  sprintf(hname3, "Entryx-Exitx_TID_%i", i + 1);
154  sprintf(hname4, "Entryx-Exitx_TEC_%i", i + 1);
155  if (pixelOutput)
156  sprintf(hname5, "Entryx-Exitx_BPIX_%i", i + 1);
157  if (pixelOutput)
158  sprintf(hname6, "Entryx-Exitx_FPIX_%i", i + 1);
159 
160  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
161  h1ex[i] = ibooker.book1D(hname1, htitle1, nbin, low[0], high[0]);
162  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
163  h2ex[i] = ibooker.book1D(hname2, htitle2, nbin, low[1], high[1]);
164  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
165  h3ex[i] = ibooker.book1D(hname3, htitle3, nbin, low[2], high[2]);
166  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
167  h4ex[i] = ibooker.book1D(hname4, htitle4, nbin, low[3], high[3]);
168  if (pixelOutput) {
169  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
170  h5ex[i] = ibooker.book1D(hname5, htitle5, nbin, low[4], high[4]);
171  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
172  h6ex[i] = ibooker.book1D(hname6, htitle6, nbin, low[5], high[5]);
173  }
174  }
175 
176  const float high0[] = {0.05, 0.06, 0.03, 0.03, 0.03, 0.03};
177  const float low0[] = {-0.05, -0.06, -0.03, -0.03, -0.03, -0.03};
178 
179  for (int i = 0; i < 12; i++) {
180  sprintf(htitle1, "Entryy-Exity in TIB %s", Region[i]);
181  sprintf(htitle2, "Entryy-Exity in TOB %s", Region[i]);
182  sprintf(htitle3, "Entryy-Exity in TID %s", Region[i]);
183  sprintf(htitle4, "Entryy-Exity in TEC %s", Region[i]);
184  if (pixelOutput)
185  sprintf(htitle5, "Entryy-Exity in BPIX %s", Region[i]);
186  if (pixelOutput)
187  sprintf(htitle6, "Entryy-Exity in FPIX %s", Region[i]);
188 
189  sprintf(hname1, "Entryy-Exity_TIB_%i", i + 1);
190  sprintf(hname2, "Entryy-Exity_TOB_%i", i + 1);
191  sprintf(hname3, "Entryy-Exity_TID_%i", i + 1);
192  sprintf(hname4, "Entryy-Exity_TEC_%i", i + 1);
193  if (pixelOutput)
194  sprintf(hname5, "Entryy-Exity_BPIX_%i", i + 1);
195  if (pixelOutput)
196  sprintf(hname6, "Entryy-Exity_FPIX_%i", i + 1);
197 
198  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
199  h1ey[i] = ibooker.book1D(hname1, htitle1, nbin, low0[0], high0[0]);
200  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
201  h2ey[i] = ibooker.book1D(hname2, htitle2, nbin, low0[1], high0[1]);
202  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
203  h3ey[i] = ibooker.book1D(hname3, htitle3, nbin, low0[2], high0[2]);
204  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
205  h4ey[i] = ibooker.book1D(hname4, htitle4, nbin, low0[3], high0[3]);
206  if (pixelOutput) {
207  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
208  h5ey[i] = ibooker.book1D(hname5, htitle5, nbin, low0[4], high0[4]);
209  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
210  h6ey[i] = ibooker.book1D(hname6, htitle6, nbin, low0[5], high0[5]);
211  }
212  }
213 
214  const float high1[] = {0.05, 0.06, 0.05, 0.06, 0.05, 0.05};
215  const float low1[] = {0., 0., 0., 0., 0., 0.};
216 
217  for (int i = 0; i < 12; i++) {
218  sprintf(htitle1, "abs(Entryz-Exitz) in TIB %s", Region[i]);
219  sprintf(htitle2, "abs(Entryz-Exitz) in TOB %s", Region[i]);
220  sprintf(htitle3, "abs(Entryz-Exitz) in TID %s", Region[i]);
221  sprintf(htitle4, "abs(Entryz-Exitz) in TEC %s", Region[i]);
222  if (pixelOutput)
223  sprintf(htitle5, "abs(Entryz-Exitz) in BPIX %s", Region[i]);
224  if (pixelOutput)
225  sprintf(htitle6, "abs(Entryz-Exitz) in FPIX %s", Region[i]);
226 
227  sprintf(hname1, "Entryz-Exitz_TIB_%i", i + 1);
228  sprintf(hname2, "Entryz-Exitz_TOB_%i", i + 1);
229  sprintf(hname3, "Entryz-Exitz_TID_%i", i + 1);
230  sprintf(hname4, "Entryz-Exitz_TEC_%i", i + 1);
231  if (pixelOutput)
232  sprintf(hname5, "Entryz-Exitz_BPIX_%i", i + 1);
233  if (pixelOutput)
234  sprintf(hname6, "Entryz-Exitz_FPIX_%i", i + 1);
235 
236  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
237  h1ez[i] = ibooker.book1D(hname1, htitle1, nbin, low1[0], high1[0]);
238  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
239  h2ez[i] = ibooker.book1D(hname2, htitle2, nbin, low1[1], high1[1]);
240  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
241  h3ez[i] = ibooker.book1D(hname3, htitle3, nbin, low1[2], high1[2]);
242  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
243  h4ez[i] = ibooker.book1D(hname4, htitle4, nbin, low1[3], high1[3]);
244  if (pixelOutput) {
245  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
246  h5ez[i] = ibooker.book1D(hname5, htitle5, nbin, low1[4], high1[4]);
247  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
248  h6ez[i] = ibooker.book1D(hname6, htitle6, nbin, low1[5], high1[5]);
249  }
250  }
251 
252  const float high2[] = {3.2, 5.0, 5.5, 6.2, 0.85, 0.5};
253  const float low2[] = {-3.2, -5.0, -5.5, -6.2, -0.85, -0.5};
254 
255  for (int i = 0; i < 12; i++) {
256  sprintf(htitle1, "Localx in TIB %s", Region[i]);
257  sprintf(htitle2, "Localx in TOB %s", Region[i]);
258  sprintf(htitle3, "Localx in TID %s", Region[i]);
259  sprintf(htitle4, "Localx in TEC %s", Region[i]);
260  if (pixelOutput)
261  sprintf(htitle5, "Localx in BPIX %s", Region[i]);
262  if (pixelOutput)
263  sprintf(htitle6, "Localx in FPIX %s", Region[i]);
264 
265  sprintf(hname1, "Localx_TIB_%i", i + 1);
266  sprintf(hname2, "Localx_TOB_%i", i + 1);
267  sprintf(hname3, "Localx_TID_%i", i + 1);
268  sprintf(hname4, "Localx_TEC_%i", i + 1);
269  if (pixelOutput)
270  sprintf(hname5, "Localx_BPIX_%i", i + 1);
271  if (pixelOutput)
272  sprintf(hname6, "Localx_FPIX_%i", i + 1);
273 
274  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
275  h1lx[i] = ibooker.book1D(hname1, htitle1, nbin, low2[0], high2[0]);
276  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
277  h2lx[i] = ibooker.book1D(hname2, htitle2, nbin, low2[1], high2[1]);
278  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
279  h3lx[i] = ibooker.book1D(hname3, htitle3, nbin, low2[2], high2[2]);
280  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
281  h4lx[i] = ibooker.book1D(hname4, htitle4, nbin, low2[3], high2[3]);
282  if (pixelOutput) {
283  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
284  h5lx[i] = ibooker.book1D(hname5, htitle5, nbin, low2[4], high2[4]);
285  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
286  h6lx[i] = ibooker.book1D(hname6, htitle6, nbin, low2[5], high2[5]);
287  }
288  }
289 
290  const float high3[] = {6.0, 10., 5.6, 10.5, 3.4, 0.52};
291  const float low3[] = {-6.0, -10., -5.6, -10.5, -3.4, -0.52};
292 
293  for (int i = 0; i < 12; i++) {
294  sprintf(htitle1, "Localy in TIB %s", Region[i]);
295  sprintf(htitle2, "Localy in TOB %s", Region[i]);
296  sprintf(htitle3, "Localy in TID %s", Region[i]);
297  sprintf(htitle4, "Localy in TEC %s", Region[i]);
298  if (pixelOutput)
299  sprintf(htitle5, "Localy in BPIX %s", Region[i]);
300  if (pixelOutput)
301  sprintf(htitle6, "Localy in FPIX %s", Region[i]);
302 
303  sprintf(hname1, "Localy_TIB_%i", i + 1);
304  sprintf(hname2, "Localy_TOB_%i", i + 1);
305  sprintf(hname3, "Localy_TID_%i", i + 1);
306  sprintf(hname4, "Localy_TEC_%i", i + 1);
307  if (pixelOutput)
308  sprintf(hname5, "Localy_BPIX_%i", i + 1);
309  if (pixelOutput)
310  sprintf(hname6, "Localy_FPIX_%i", i + 1);
311 
312  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIBHit");
313  h1ly[i] = ibooker.book1D(hname1, htitle1, nbin, low3[0], high3[0]);
314  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TOBHit");
315  h2ly[i] = ibooker.book1D(hname2, htitle2, nbin, low3[1], high3[1]);
316  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TIDHit");
317  h3ly[i] = ibooker.book1D(hname3, htitle3, nbin, low3[2], high3[2]);
318  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/TECHit");
319  h4ly[i] = ibooker.book1D(hname4, htitle4, nbin, low3[3], high3[3]);
320  if (pixelOutput) {
321  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/BPIXHit");
322  h5ly[i] = ibooker.book1D(hname5, htitle5, nbin, low3[4], high3[4]);
323  ibooker.setCurrentFolder("TrackerHitsV/TrackerHit/FPIXHit");
324  h6ly[i] = ibooker.book1D(hname6, htitle6, nbin, low3[5], high3[5]);
325  }
326  }
327  }
328 }
329 
331  // don't try to delete any pointers - they're handled by DQM machinery
332 }
333 
335  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
336 
337  // iterator to access containers
338  edm::PSimHitContainer::const_iterator itHit;
340  // get Pixel Barrel information
342  // extract low container
343  edm::Handle<edm::PSimHitContainer> PxlBrlLowContainer;
344  e.getByToken(edmPSimHitContainer_pxlBrlLow_Token_, PxlBrlLowContainer);
345  if (!PxlBrlLowContainer.isValid()) {
346  edm::LogError("TrackerHitAnalyzer::analyze") << "Unable to find TrackerHitsPixelBarrelLowTof in event!";
347  return;
348  }
349  // extract high container
350  edm::Handle<edm::PSimHitContainer> PxlBrlHighContainer;
351  e.getByToken(edmPSimHitContainer_pxlBrlHigh_Token_, PxlBrlHighContainer);
352  if (!PxlBrlHighContainer.isValid()) {
353  edm::LogError("TrackerHitAnalyzer::analyze") << "Unable to find TrackerHitsPixelBarrelHighTof in event!";
354  return;
355  }
357  // get Pixel Forward information
359  // extract low container
360  edm::Handle<edm::PSimHitContainer> PxlFwdLowContainer;
361  e.getByToken(edmPSimHitContainer_pxlFwdLow_Token_, PxlFwdLowContainer);
362  if (!PxlFwdLowContainer.isValid()) {
363  edm::LogError("TrackerHitAnalyzer::analyze") << "Unable to find TrackerHitsPixelEndcapLowTof in event!";
364  return;
365  }
366  // extract high container
367  edm::Handle<edm::PSimHitContainer> PxlFwdHighContainer;
368  e.getByToken(edmPSimHitContainer_pxlFwdHigh_Token_, PxlFwdHighContainer);
369  if (!PxlFwdHighContainer.isValid()) {
370  edm::LogError("TrackerHitAnalyzer::analyze") << "Unable to find TrackerHitsPixelEndcapHighTof in event!";
371  return;
372  }
374  // get Silicon TIB information
376  // extract TIB low container
377  edm::Handle<edm::PSimHitContainer> SiTIBLowContainer;
378  e.getByToken(edmPSimHitContainer_siTIBLow_Token_, SiTIBLowContainer);
379  if (!SiTIBLowContainer.isValid()) {
380  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTIBLowTof in event!";
381  return;
382  }
384  // extract TIB high container
385  edm::Handle<edm::PSimHitContainer> SiTIBHighContainer;
386  e.getByToken(edmPSimHitContainer_siTIBHigh_Token_, SiTIBHighContainer);
387  if (!SiTIBHighContainer.isValid()) {
388  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTIBHighTof in event!";
389  return;
390  }
392  // get Silicon TOB information
394  // extract TOB low container
395  edm::Handle<edm::PSimHitContainer> SiTOBLowContainer;
396  e.getByToken(edmPSimHitContainer_siTOBLow_Token_, SiTOBLowContainer);
397  if (!SiTOBLowContainer.isValid()) {
398  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTOBLowTof in event!";
399  return;
400  }
402  // extract TOB high container
403  edm::Handle<edm::PSimHitContainer> SiTOBHighContainer;
404  e.getByToken(edmPSimHitContainer_siTOBHigh_Token_, SiTOBHighContainer);
405  if (!SiTOBHighContainer.isValid()) {
406  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTOBHighTof in event!";
407  return;
408  }
409 
411  // get Silicon TID information
413  // extract TID low container
414  edm::Handle<edm::PSimHitContainer> SiTIDLowContainer;
415  e.getByToken(edmPSimHitContainer_siTIDLow_Token_, SiTIDLowContainer);
416  if (!SiTIDLowContainer.isValid()) {
417  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTIDLowTof in event!";
418  return;
419  }
421  // extract TID high container
422  edm::Handle<edm::PSimHitContainer> SiTIDHighContainer;
423  e.getByToken(edmPSimHitContainer_siTIDHigh_Token_, SiTIDHighContainer);
424  if (!SiTIDHighContainer.isValid()) {
425  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTIDHighTof in event!";
426  return;
427  }
429  // get Silicon TEC information
431  // extract TEC low container
432  edm::Handle<edm::PSimHitContainer> SiTECLowContainer;
433  e.getByToken(edmPSimHitContainer_siTECLow_Token_, SiTECLowContainer);
434  if (!SiTECLowContainer.isValid()) {
435  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTECLowTof in event!";
436  return;
437  }
439  // extract TEC high container
440  edm::Handle<edm::PSimHitContainer> SiTECHighContainer;
441  e.getByToken(edmPSimHitContainer_siTECHigh_Token_, SiTECHighContainer);
442  if (!SiTECHighContainer.isValid()) {
443  edm::LogError("TrackerHitProducer::analyze") << "Unable to find TrackerHitsTECHighTof in event!";
444  return;
445  }
446 
448  // get G4Track information
450 
452  e.getByToken(edmSimTrackContainerToken_, G4TrkContainer);
453  if (!G4TrkContainer.isValid()) {
454  edm::LogError("TrackerHitAnalyzer::analyze") << "Unable to find SimTrack in event!";
455  return;
456  }
457 
458  // Get geometry information
459  const auto &tracker = &c.getData(tGeomEsToken_);
460 
461  int ir = -100;
462  edm::SimTrackContainer::const_iterator itTrk;
463  for (itTrk = G4TrkContainer->begin(); itTrk != G4TrkContainer->end(); ++itTrk) {
464  // std::cout << "itTrk = "<< itTrk << std::endl;
465  double eta = 0, p = 0;
466  const CLHEP::HepLorentzVector &G4Trk = CLHEP::HepLorentzVector(
467  itTrk->momentum().x(), itTrk->momentum().y(), itTrk->momentum().z(), itTrk->momentum().e());
468  p = sqrt(G4Trk[0] * G4Trk[0] + G4Trk[1] * G4Trk[1] + G4Trk[2] * G4Trk[2]);
469  if (p == 0)
470  edm::LogError("TrackerHitAnalyzer::analyze") << "TrackerTest::INFO: Primary has p = 0 ";
471  else {
472  double costheta = G4Trk[2] / p;
473  double theta = acos(TMath::Min(TMath::Max(costheta, -1.), 1.));
474  eta = -log(tan(theta / 2));
475 
476  if (eta > 0.0 && eta <= 0.5)
477  ir = 0;
478  if (eta > 0.5 && eta <= 1.0)
479  ir = 1;
480  if (eta > 1.0 && eta <= 1.5)
481  ir = 2;
482  if (eta > 1.5 && eta <= 2.0)
483  ir = 3;
484  if (eta > 2.0 && eta <= 2.5)
485  ir = 4;
486  if (eta > 2.5)
487  ir = 5;
488 
489  if (eta > -0.5 && eta <= 0.0)
490  ir = 6;
491  if (eta > -1.0 && eta <= -0.5)
492  ir = 7;
493  if (eta > -1.5 && eta <= -1.0)
494  ir = 8;
495  if (eta > -2.0 && eta <= -1.5)
496  ir = 9;
497  if (eta > -2.5 && eta <= -2.0)
498  ir = 10;
499  if (eta <= -2.5)
500  ir = 11;
501  // edm::LogInfo("EventInfo") << " eta = " << eta << " ir = " <<
502  // ir;
503  // std::cout << " " << std::endl;
504  // std::cout << "eta " << eta << " ir = " << ir << std::endl;
505  // std::cout << " " << std::endl;
506  }
507  }
509  // get Pixel information
511  // If Phase 1, do not run - will run in new Phase 1 module
512 
513  if (pixelOutput) {
514  for (itHit = PxlBrlLowContainer->begin(); itHit != PxlBrlLowContainer->end(); ++itHit) {
515  DetId detid = DetId(itHit->detUnitId());
516  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
517  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
518  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
519  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
520  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
521  htofz->Fill(gpos.z(), itHit->timeOfFlight());
522 
523  h5e[ir]->Fill(itHit->energyLoss());
524  h5ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
525  h5ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
526  h5ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
527  h5lx[ir]->Fill(itHit->localPosition().x());
528  h5ly[ir]->Fill(itHit->localPosition().y());
529  }
530  for (itHit = PxlBrlHighContainer->begin(); itHit != PxlBrlHighContainer->end(); ++itHit) {
531  DetId detid = DetId(itHit->detUnitId());
532  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
533  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
534  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
535  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
536  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
537  htofz->Fill(gpos.z(), itHit->timeOfFlight());
538 
539  h5e[ir]->Fill(itHit->energyLoss());
540  h5ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
541  h5ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
542  h5ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
543  h5lx[ir]->Fill(itHit->localPosition().x());
544  h5ly[ir]->Fill(itHit->localPosition().y());
545  }
546  for (itHit = PxlFwdLowContainer->begin(); itHit != PxlFwdLowContainer->end(); ++itHit) {
547  DetId detid = DetId(itHit->detUnitId());
548  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
549  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
550  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
551  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
552  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
553  htofz->Fill(gpos.z(), itHit->timeOfFlight());
554 
555  h6e[ir]->Fill(itHit->energyLoss());
556  h6ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
557  h6ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
558  h6ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
559  h6lx[ir]->Fill(itHit->localPosition().x());
560  h6ly[ir]->Fill(itHit->localPosition().y());
561  }
562  for (itHit = PxlFwdHighContainer->begin(); itHit != PxlFwdHighContainer->end(); ++itHit) {
563  DetId detid = DetId(itHit->detUnitId());
564  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
565  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
566  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
567  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
568  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
569  htofz->Fill(gpos.z(), itHit->timeOfFlight());
570 
571  h6e[ir]->Fill(itHit->energyLoss());
572  h6ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
573  h6ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
574  h6ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
575  h6lx[ir]->Fill(itHit->localPosition().x());
576  h6ly[ir]->Fill(itHit->localPosition().y());
577  }
578  }
580  // get TIB information
582  for (itHit = SiTIBLowContainer->begin(); itHit != SiTIBLowContainer->end(); ++itHit) {
583  DetId detid = DetId(itHit->detUnitId());
584  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
585  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
586  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
587  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
588  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
589  htofz->Fill(gpos.z(), itHit->timeOfFlight());
590 
591  h1e[ir]->Fill(itHit->energyLoss());
592  h1ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
593  h1ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
594  h1ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
595  h1lx[ir]->Fill(itHit->localPosition().x());
596  h1ly[ir]->Fill(itHit->localPosition().y());
597  }
598  for (itHit = SiTIBHighContainer->begin(); itHit != SiTIBHighContainer->end(); ++itHit) {
599  DetId detid = DetId(itHit->detUnitId());
600  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
601  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
602  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
603  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
604  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
605  htofz->Fill(gpos.z(), itHit->timeOfFlight());
606 
607  h1e[ir]->Fill(itHit->energyLoss());
608  h1ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
609  h1ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
610  h1ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
611  h1lx[ir]->Fill(itHit->localPosition().x());
612  h1ly[ir]->Fill(itHit->localPosition().y());
613  }
615  // get TOB information
617  for (itHit = SiTOBLowContainer->begin(); itHit != SiTOBLowContainer->end(); ++itHit) {
618  DetId detid = DetId(itHit->detUnitId());
619  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
620  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
621  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
622  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
623  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
624  htofz->Fill(gpos.z(), itHit->timeOfFlight());
625 
626  h2e[ir]->Fill(itHit->energyLoss());
627  h2ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
628  h2ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
629  h2ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
630  h2lx[ir]->Fill(itHit->localPosition().x());
631  h2ly[ir]->Fill(itHit->localPosition().y());
632  }
633  for (itHit = SiTOBHighContainer->begin(); itHit != SiTOBHighContainer->end(); ++itHit) {
634  DetId detid = DetId(itHit->detUnitId());
635  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
636  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
637  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
638  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
639  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
640  htofz->Fill(gpos.z(), itHit->timeOfFlight());
641 
642  h2e[ir]->Fill(itHit->energyLoss());
643  h2ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
644  h2ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
645  h2ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
646  h2lx[ir]->Fill(itHit->localPosition().x());
647  h2ly[ir]->Fill(itHit->localPosition().y());
648  }
650  // get TID information
652  for (itHit = SiTIDLowContainer->begin(); itHit != SiTIDLowContainer->end(); ++itHit) {
653  DetId detid = DetId(itHit->detUnitId());
654  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
655  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
656  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
657  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
658  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
659  htofz->Fill(gpos.z(), itHit->timeOfFlight());
660 
661  h3e[ir]->Fill(itHit->energyLoss());
662  h3ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
663  h3ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
664  h3ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
665  h3lx[ir]->Fill(itHit->localPosition().x());
666  h3ly[ir]->Fill(itHit->localPosition().y());
667  }
668  for (itHit = SiTIDHighContainer->begin(); itHit != SiTIDHighContainer->end(); ++itHit) {
669  DetId detid = DetId(itHit->detUnitId());
670  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
671  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
672  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
673  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
674  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
675  htofz->Fill(gpos.z(), itHit->timeOfFlight());
676 
677  h3e[ir]->Fill(itHit->energyLoss());
678  h3ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
679  h3ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
680  h3ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
681  h3lx[ir]->Fill(itHit->localPosition().x());
682  h3ly[ir]->Fill(itHit->localPosition().y());
683  }
685  // get TEC information
687  for (itHit = SiTECLowContainer->begin(); itHit != SiTECLowContainer->end(); ++itHit) {
688  DetId detid = DetId(itHit->detUnitId());
689  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
690  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
691  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
692  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
693  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
694  htofz->Fill(gpos.z(), itHit->timeOfFlight());
695 
696  h4e[ir]->Fill(itHit->energyLoss());
697  h4ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
698  h4ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
699  h4ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
700  h4lx[ir]->Fill(itHit->localPosition().x());
701  h4ly[ir]->Fill(itHit->localPosition().y());
702  }
703  for (itHit = SiTECHighContainer->begin(); itHit != SiTECHighContainer->end(); ++itHit) {
704  DetId detid = DetId(itHit->detUnitId());
705  const GeomDetUnit *det = (const GeomDetUnit *)tracker->idToDetUnit(detid);
706  GlobalPoint gpos = det->toGlobal(itHit->localPosition());
707  htofeta->Fill(gpos.eta(), itHit->timeOfFlight());
708  htofphi->Fill(gpos.phi().degrees(), itHit->timeOfFlight());
709  htofr->Fill(gpos.mag(), itHit->timeOfFlight());
710  htofz->Fill(gpos.z(), itHit->timeOfFlight());
711 
712  h4e[ir]->Fill(itHit->energyLoss());
713  h4ex[ir]->Fill(itHit->entryPoint().x() - itHit->exitPoint().x());
714  h4ey[ir]->Fill(itHit->entryPoint().y() - itHit->exitPoint().y());
715  h4ez[ir]->Fill(std::fabs(itHit->entryPoint().z() - itHit->exitPoint().z()));
716  h4lx[ir]->Fill(itHit->localPosition().x());
717  h4ly[ir]->Fill(itHit->localPosition().y());
718  }
719 
720  return;
721 }
722 
MonitorElement * h1e[12]
void bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &run, const edm::EventSetup &es) override
MonitorElement * h2lx[12]
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
MonitorElement * h4e[12]
MonitorElement * h1lx[12]
MonitorElement * h6ex[12]
MonitorElement * h1ey[12]
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTECLow_Token_
MonitorElement * h3ly[12]
MonitorElement * h5ez[12]
MonitorElement * h4ex[12]
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
edm::EDGetTokenT< edm::SimTrackContainer > edmSimTrackContainerToken_
T z() const
Definition: PV3DBase.h:61
MonitorElement * h1ez[12]
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
~TrackerHitAnalyzer() override
Destructor.
MonitorElement * h3ez[12]
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_pxlBrlLow_Token_
T eta() const
Definition: PV3DBase.h:73
MonitorElement * h5ey[12]
MonitorElement * h3lx[12]
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTOBLow_Token_
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_pxlFwdLow_Token_
Log< level::Error, false > LogError
MonitorElement * htofr
MonitorElement * h6ly[12]
MonitorElement * h5ex[12]
MonitorElement * h5ly[12]
MonitorElement * htofeta
void Fill(long long x)
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTIBLow_Token_
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Analyze.
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTOBHigh_Token_
T sqrt(T t)
Definition: SSEVec.h:19
MonitorElement * h1ex[12]
MonitorElement * h2e[12]
T mag() const
Definition: PV3DBase.h:64
MonitorElement * h6ey[12]
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
MonitorElement * h4lx[12]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * h2ex[12]
MonitorElement * h6ez[12]
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTIDLow_Token_
MonitorElement * h2ez[12]
MonitorElement * htofz
Log< level::Info, false > LogInfo
Definition: DetId.h:17
MonitorElement * htofphi
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_pxlBrlHigh_Token_
MonitorElement * h5lx[12]
MonitorElement * h4ey[12]
MonitorElement * h1ly[12]
T1 degrees() const
Definition: Phi.h:107
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomEsToken_
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
MonitorElement * h4ly[12]
bool isValid() const
Definition: HandleBase.h:70
MonitorElement * h3ex[12]
MonitorElement * h3ey[12]
HLT enums.
MonitorElement * h2ey[12]
MonitorElement * h2ly[12]
MonitorElement * h3e[12]
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTIBHigh_Token_
std::vector< PSimHit > PSimHitContainer
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_pxlFwdHigh_Token_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
TrackerHitAnalyzer(const edm::ParameterSet &ps)
Constructor.
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTECHigh_Token_
MonitorElement * h6e[12]
Geom::Theta< T > theta() const
std::vector< SimTrack > SimTrackContainer
edm::EDGetTokenT< edm::PSimHitContainer > edmPSimHitContainer_siTIDHigh_Token_
MonitorElement * h5e[12]
MonitorElement * h6lx[12]
Definition: Run.h:45
MonitorElement * h4ez[12]