CMS 3D CMS Logo

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