CMS 3D CMS Logo

SiPixelRecHitsValid.cc
Go to the documentation of this file.
1 // SiPixelRecHitsValid.cc
2 // Description: see SiPixelRecHitsValid.h
3 // Author: Jason Shaev, JHU
4 // Created 6/7/06
5 //
6 // G. Giurgiu, JHU (ggiurgiu@pha.jhu.edu)
7 // added pull distributions (12/27/06)
8 //--------------------------------
9 
11 
23 
25 
32 
42 
45 
46 #include <cmath>
47 
49  : trackerHitAssociatorConfig_(ps, consumesCollector()),
50  siPixelRecHitCollectionToken_(consumes<SiPixelRecHitCollection>(ps.getParameter<edm::InputTag>("src"))) {}
51 
53 
55 
57  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustBPIX");
58 
59  Char_t histo[200];
60 
61  // ---------------------------------------------------------------
62  // All histograms that depend on plaquette number have 7 indexes.
63  // The first 4 (0-3) correspond to Panel 1 plaquettes 1-4.
64  // The last 3 (4-6) correspond to Panel 2 plaquettes 1-3.
65  // ---------------------------------------------------------------
66 
67  //Cluster y-size by module number for barrel
68  for (int i = 0; i < 8; i++) {
69  sprintf(histo, "Clust_y_size_Module%d", i + 1);
70  clustYSizeModule[i] = ibooker.book1D(histo, "Cluster y-size by Module", 20, 0.5, 20.5);
71  } // end for
72 
73  //Cluster x-size by layer for barrel
74  for (int i = 0; i < 3; i++) {
75  sprintf(histo, "Clust_x_size_Layer%d", i + 1);
76  clustXSizeLayer[i] = ibooker.book1D(histo, "Cluster x-size by Layer", 20, 0.5, 20.5);
77  } // end for
78 
79  //Cluster charge by module for 3 layers of barrel
80  for (int i = 0; i < 8; i++) {
81  //Cluster charge by module for Layer1
82  sprintf(histo, "Clust_charge_Layer1_Module%d", i + 1);
83  clustChargeLayer1Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 1 by Module", 50, 0., 200000.);
84 
85  //Cluster charge by module for Layer2
86  sprintf(histo, "Clust_charge_Layer2_Module%d", i + 1);
87  clustChargeLayer2Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 2 by Module", 50, 0., 200000.);
88 
89  //Cluster charge by module for Layer3
90  sprintf(histo, "Clust_charge_Layer3_Module%d", i + 1);
91  clustChargeLayer3Modules[i] = ibooker.book1D(histo, "Cluster charge Layer 3 by Module", 50, 0., 200000.);
92  } // end for
93 
94  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/clustFPIX");
95  //Cluster x-size, y-size, and charge by plaquette for Disks in Forward
96  for (int i = 0; i < 7; i++) {
97  //Cluster x-size for Disk1 by Plaquette
98  sprintf(histo, "Clust_x_size_Disk1_Plaquette%d", i + 1);
99  clustXSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk1 by Plaquette", 20, 0.5, 20.5);
100 
101  //Cluster x-size for Disk2 by Plaquette
102  sprintf(histo, "Clust_x_size_Disk2_Plaquette%d", i + 1);
103  clustXSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster X-size for Disk2 by Plaquette", 20, 0.5, 20.5);
104 
105  //Cluster y-size for Disk1 by Plaquette
106  sprintf(histo, "Clust_y_size_Disk1_Plaquette%d", i + 1);
107  clustYSizeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk1 by Plaquette", 20, 0.5, 20.5);
108 
109  //Cluster y-size for Disk2 by Plaquette
110  sprintf(histo, "Clust_y_size_Disk2_Plaquette%d", i + 1);
111  clustYSizeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster Y-size for Disk2 by Plaquette", 20, 0.5, 20.5);
112 
113  //Cluster charge for Disk1 by Plaquette
114  sprintf(histo, "Clust_charge_Disk1_Plaquette%d", i + 1);
115  clustChargeDisk1Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk1 by Plaquette", 50, 0., 200000.);
116 
117  //Cluster charge for Disk2 by Plaquette
118  sprintf(histo, "Clust_charge_Disk2_Plaquette%d", i + 1);
119  clustChargeDisk2Plaquettes[i] = ibooker.book1D(histo, "Cluster charge for Disk2 by Plaquette", 50, 0., 200000.);
120  } // end for
121 
122  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitBPIX");
123  //RecHit Bunch crossing all barrel hits
124  recHitBunchB = ibooker.book1D("RecHit_Bunch_Barrel", "RecHit Bunch Crossing, Barrel", 20, -10., 10.);
125 
126  //RecHit Event, in-time bunch, all barrel hits
127  recHitEventB = ibooker.book1D("RecHit_Event_Barrel", "RecHit Event (in-time bunch), Barrel", 100, 0., 100.);
128 
129  //RecHit X Resolution all barrel hits
130  recHitXResAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Res All Modules in Barrel", 100, -200., 200.);
131 
132  //RecHit Y Resolution all barrel hits
133  recHitYResAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Res All Modules in Barrel", 100, -200., 200.);
134 
135  //RecHit X distribution for full modules for barrel
136  recHitXFullModules = ibooker.book1D("RecHit_x_FullModules", "RecHit X distribution for full modules", 100, -2., 2.);
137 
138  //RecHit X distribution for half modules for barrel
139  recHitXHalfModules = ibooker.book1D("RecHit_x_HalfModules", "RecHit X distribution for half modules", 100, -1., 1.);
140 
141  //RecHit Y distribution all modules for barrel
142  recHitYAllModules = ibooker.book1D("RecHit_y_AllModules", "RecHit Y distribution for all modules", 100, -4., 4.);
143 
144  //RecHit X resolution for flipped and unflipped ladders by layer for barrel
145  for (int i = 0; i < 3; i++) {
146  //RecHit no. of matched simHits all ladders by layer
147  sprintf(histo, "RecHit_NsimHit_Layer%d", i + 1);
148  recHitNsimHitLayer[i] = ibooker.book1D(histo, "RecHit Number of simHits by Layer", 30, 0., 30.);
149 
150  //RecHit X resolution for flipped ladders by layer
151  sprintf(histo, "RecHit_XRes_FlippedLadder_Layer%d", i + 1);
152  recHitXResFlippedLadderLayers[i] = ibooker.book1D(histo, "RecHit XRes Flipped Ladders by Layer", 100, -200., 200.);
153 
154  //RecHit X resolution for unflipped ladders by layer
155  sprintf(histo, "RecHit_XRes_UnFlippedLadder_Layer%d", i + 1);
157  ibooker.book1D(histo, "RecHit XRes NonFlipped Ladders by Layer", 100, -200., 200.);
158  } // end for
159 
160  //RecHit Y resolutions for layers by module for barrel
161  for (int i = 0; i < 8; i++) {
162  //Rec Hit Y resolution by module for Layer1
163  sprintf(histo, "RecHit_YRes_Layer1_Module%d", i + 1);
164  recHitYResLayer1Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer1 by module", 100, -200., 200.);
165 
166  //RecHit Y resolution by module for Layer2
167  sprintf(histo, "RecHit_YRes_Layer2_Module%d", i + 1);
168  recHitYResLayer2Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer2 by module", 100, -200., 200.);
169 
170  //RecHit Y resolution by module for Layer3
171  sprintf(histo, "RecHit_YRes_Layer3_Module%d", i + 1);
172  recHitYResLayer3Modules[i] = ibooker.book1D(histo, "RecHit YRes Layer3 by module", 100, -200., 200.);
173  } // end for
174 
175  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitFPIX");
176  //RecHit Bunch crossing all plaquettes
177  recHitBunchF = ibooker.book1D("RecHit_Bunch_Forward", "RecHit Bunch Crossing, Forward", 20, -10., 10.);
178 
179  //RecHit Event, in-time bunch, all plaquettes
180  recHitEventF = ibooker.book1D("RecHit_Event_Forward", "RecHit Event (in-time bunch), Forward", 100, 0., 100.);
181 
182  //RecHit No. of simHits, by disk
183  recHitNsimHitDisk1 = ibooker.book1D("RecHit_NsimHit_Disk1", "RecHit Number of simHits, Disk1", 30, 0., 30.);
184  recHitNsimHitDisk2 = ibooker.book1D("RecHit_NsimHit_Disk2", "RecHit Number of simHits, Disk2", 30, 0., 30.);
185 
186  //RecHit X resolution all plaquettes
187  recHitXResAllF = ibooker.book1D("RecHit_xres_f_All", "RecHit X Res All in Forward", 100, -200., 200.);
188 
189  //RecHit Y resolution all plaquettes
190  recHitYResAllF = ibooker.book1D("RecHit_yres_f_All", "RecHit Y Res All in Forward", 100, -200., 200.);
191 
192  //RecHit X distribution for plaquette with x-size 1 in forward
194  ibooker.book1D("RecHit_x_Plaquette_xsize1", "RecHit X Distribution for plaquette x-size1", 100, -2., 2.);
195 
196  //RecHit X distribution for plaquette with x-size 2 in forward
198  ibooker.book1D("RecHit_x_Plaquette_xsize2", "RecHit X Distribution for plaquette x-size2", 100, -2., 2.);
199 
200  //RecHit Y distribution for plaquette with y-size 2 in forward
202  ibooker.book1D("RecHit_y_Plaquette_ysize2", "RecHit Y Distribution for plaquette y-size2", 100, -4., 4.);
203 
204  //RecHit Y distribution for plaquette with y-size 3 in forward
206  ibooker.book1D("RecHit_y_Plaquette_ysize3", "RecHit Y Distribution for plaquette y-size3", 100, -4., 4.);
207 
208  //RecHit Y distribution for plaquette with y-size 4 in forward
210  ibooker.book1D("RecHit_y_Plaquette_ysize4", "RecHit Y Distribution for plaquette y-size4", 100, -4., 4.);
211 
212  //RecHit Y distribution for plaquette with y-size 5 in forward
214  ibooker.book1D("RecHit_y_Plaquette_ysize5", "RecHit Y Distribution for plaquette y-size5", 100, -4., 4.);
215 
216  //X and Y resolutions for both disks by plaquette in forward
217  for (int i = 0; i < 7; i++) {
218  //X resolution for Disk1 by plaquette
219  sprintf(histo, "RecHit_XRes_Disk1_Plaquette%d", i + 1);
220  recHitXResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk1 by plaquette", 100, -200., 200.);
221  //X resolution for Disk2 by plaquette
222  sprintf(histo, "RecHit_XRes_Disk2_Plaquette%d", i + 1);
223  recHitXResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XRes Disk2 by plaquette", 100, -200., 200.);
224 
225  //Y resolution for Disk1 by plaquette
226  sprintf(histo, "RecHit_YRes_Disk1_Plaquette%d", i + 1);
227  recHitYResDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk1 by plaquette", 100, -200., 200.);
228  //Y resolution for Disk2 by plaquette
229  sprintf(histo, "RecHit_YRes_Disk2_Plaquette%d", i + 1);
230  recHitYResDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YRes Disk2 by plaquette", 100, -200., 200.);
231  }
232 
233  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsBPIX");
234  recHitXPullAllB = ibooker.book1D("RecHit_xres_b_All", "RecHit X Pull All Modules in Barrel", 100, -10.0, 10.0);
235  recHitYPullAllB = ibooker.book1D("RecHit_yres_b_All", "RecHit Y Pull All Modules in Barrel", 100, -10.0, 10.0);
236 
237  for (int i = 0; i < 3; i++) {
238  sprintf(histo, "RecHit_XPull_FlippedLadder_Layer%d", i + 1);
240  ibooker.book1D(histo, "RecHit XPull Flipped Ladders by Layer", 100, -10.0, 10.0);
241 
242  sprintf(histo, "RecHit_XPull_UnFlippedLadder_Layer%d", i + 1);
244  ibooker.book1D(histo, "RecHit XPull NonFlipped Ladders by Layer", 100, -10.0, 10.0);
245  }
246 
247  for (int i = 0; i < 8; i++) {
248  sprintf(histo, "RecHit_YPull_Layer1_Module%d", i + 1);
249  recHitYPullLayer1Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer1 by module", 100, -10.0, 10.0);
250 
251  sprintf(histo, "RecHit_YPull_Layer2_Module%d", i + 1);
252  recHitYPullLayer2Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer2 by module", 100, -10.0, 10.0);
253 
254  sprintf(histo, "RecHit_YPull_Layer3_Module%d", i + 1);
255  recHitYPullLayer3Modules[i] = ibooker.book1D(histo, "RecHit YPull Layer3 by module", 100, -10.0, 10.0);
256  }
257 
258  ibooker.setCurrentFolder("TrackerRecHitsV/TrackerRecHits/Pixel/recHitPullsFPIX");
259  recHitXPullAllF = ibooker.book1D("RecHit_XPull_f_All", "RecHit X Pull All in Forward", 100, -10.0, 10.0);
260 
261  recHitYPullAllF = ibooker.book1D("RecHit_YPull_f_All", "RecHit Y Pull All in Forward", 100, -10.0, 10.0);
262 
263  for (int i = 0; i < 7; i++) {
264  sprintf(histo, "RecHit_XPull_Disk1_Plaquette%d", i + 1);
265  recHitXPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk1 by plaquette", 100, -10.0, 10.0);
266  sprintf(histo, "RecHit_XPull_Disk2_Plaquette%d", i + 1);
267  recHitXPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit XPull Disk2 by plaquette", 100, -10.0, 10.0);
268 
269  sprintf(histo, "RecHit_YPull_Disk1_Plaquette%d", i + 1);
270  recHitYPullDisk1Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk1 by plaquette", 100, -10.0, 10.0);
271 
272  sprintf(histo, "RecHit_YPull_Disk2_Plaquette%d", i + 1);
273  recHitYPullDisk2Plaquettes[i] = ibooker.book1D(histo, "RecHit YPull Disk2 by plaquette", 100, -10.0, 10.0);
274  }
275 }
276 
278  //Retrieve tracker topology from geometry
280  es.get<TrackerTopologyRcd>().get(tTopoHand);
281  const TrackerTopology* tTopo = tTopoHand.product();
282 
283  edm::LogInfo("EventInfo") << " Run = " << e.id().run() << " Event = " << e.id().event();
284  if (e.id().event() % 1000 == 0)
285  std::cout << " Run = " << e.id().run() << " Event = " << e.id().event() << std::endl;
286 
287  //Get RecHits
290 
291  //Get event setup
293  es.get<TrackerDigiGeometryRecord>().get(geom);
294  const TrackerGeometry& theTracker(*geom);
295 
297 
298  //iterate over detunits
299  for (TrackerGeometry::DetContainer::const_iterator it = geom->dets().begin(); it != geom->dets().end(); it++) {
300  DetId detId = ((*it)->geographicalId());
301  unsigned int subid = detId.subdetId();
302 
303  if (!((subid == 1) || (subid == 2)))
304  continue;
305 
306  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>(theTracker.idToDet(detId));
307 
308  SiPixelRecHitCollection::const_iterator pixeldet = recHitColl->find(detId);
309  if (pixeldet == recHitColl->end())
310  continue;
311  SiPixelRecHitCollection::DetSet pixelrechitRange = *pixeldet;
312  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
313  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.end();
314  SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
315  std::vector<PSimHit> matched;
316 
317  //----Loop over rechits for this detId
318  for (; pixeliter != pixelrechitRangeIteratorEnd; pixeliter++) {
319  matched.clear();
320  matched = associate.associateHit(*pixeliter);
321 
322  if (!matched.empty()) {
323  float closest = 9999.9;
324  std::vector<PSimHit>::const_iterator closestit = matched.begin();
325  LocalPoint lp = pixeliter->localPosition();
326  float rechit_x = lp.x();
327  float rechit_y = lp.y();
328 
329  //loop over sim hits and fill closet
330  for (std::vector<PSimHit>::const_iterator m = matched.begin(); m < matched.end(); m++) {
331  float sim_x1 = (*m).entryPoint().x();
332  float sim_x2 = (*m).exitPoint().x();
333  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
334 
335  float sim_y1 = (*m).entryPoint().y();
336  float sim_y2 = (*m).exitPoint().y();
337  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
338 
339  float x_res = fabs(sim_xpos - rechit_x);
340  float y_res = fabs(sim_ypos - rechit_y);
341 
342  float dist = sqrt(x_res * x_res + y_res * y_res);
343 
344  if (dist < closest) {
345  closest = dist;
346  closestit = m;
347  }
348  } // end sim hit loop
349 
350  if (subid == 1) { //<----------barrel
351  fillBarrel(*pixeliter, *closestit, detId, theGeomDet, tTopo);
352  } // end barrel
353  if (subid == 2) { // <-------forward
354  fillForward(*pixeliter, *closestit, detId, theGeomDet, tTopo);
355  }
356 
357  } // end matched emtpy
358 
359  int NsimHit = matched.size();
360  if (subid == 1) { //<----------barrel
361  for (unsigned int i = 0; i < 3; i++)
362  if (tTopo->pxbLayer(detId) == i + 1)
363  recHitNsimHitLayer[i]->Fill(NsimHit);
364  } // end barrel
365  if (subid == 2) { // <-------forward
366  if (tTopo->pxfDisk(detId) == 1)
367  recHitNsimHitDisk1->Fill(NsimHit);
368  else
369  recHitNsimHitDisk2->Fill(NsimHit);
370  }
371  } // <-----end rechit loop
372  } // <------ end detunit loop
373 }
374 
376  const PSimHit& simHit,
377  DetId detId,
378  const PixelGeomDetUnit* theGeomDet,
379  const TrackerTopology* tTopo) {
380  const float cmtomicron = 10000.0;
381 
382  int bunch = simHit.eventId().bunchCrossing();
383  int event = simHit.eventId().event();
384 
385  recHitBunchB->Fill(bunch);
386  if (bunch == 0)
388 
389  LocalPoint lp = recHit.localPosition();
390  float lp_y = lp.y();
391  float lp_x = lp.x();
392 
393  LocalError lerr = recHit.localPositionError();
394  float lerr_x = sqrt(lerr.xx());
395  float lerr_y = sqrt(lerr.yy());
396 
397  recHitYAllModules->Fill(lp_y);
398 
399  float sim_x1 = simHit.entryPoint().x();
400  float sim_x2 = simHit.exitPoint().x();
401  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
402  float res_x = (lp.x() - sim_xpos) * cmtomicron;
403 
404  recHitXResAllB->Fill(res_x);
405 
406  float sim_y1 = simHit.entryPoint().y();
407  float sim_y2 = simHit.exitPoint().y();
408  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
409  float res_y = (lp.y() - sim_ypos) * cmtomicron;
410 
411  recHitYResAllB->Fill(res_y);
412 
413  float pull_x = (lp_x - sim_xpos) / lerr_x;
414  float pull_y = (lp_y - sim_ypos) / lerr_y;
415 
416  recHitXPullAllB->Fill(pull_x);
417  recHitYPullAllB->Fill(pull_y);
418 
419  int rows = theGeomDet->specificTopology().nrows();
420 
421  if (rows == 160) {
422  recHitXFullModules->Fill(lp_x);
423  } else if (rows == 80) {
424  recHitXHalfModules->Fill(lp_x);
425  }
426 
427  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 0.)).perp();
428  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0., 0., 1.)).perp();
429 
430  if (tmp2 < tmp1) { // flipped
431  for (unsigned int i = 0; i < 3; i++) {
432  if (tTopo->pxbLayer(detId) == i + 1) {
435  }
436  }
437  } else {
438  for (unsigned int i = 0; i < 3; i++) {
439  if (tTopo->pxbLayer(detId) == i + 1) {
442  }
443  }
444  }
445 
446  //get cluster
447  SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
448 
449  // fill module dependent info
450  for (unsigned int i = 0; i < 8; i++) {
451  if (tTopo->pxbModule(detId) == i + 1) {
452  int sizeY = (*clust).sizeY();
453  clustYSizeModule[i]->Fill(sizeY);
454 
455  if (tTopo->pxbLayer(detId) == 1) {
456  float charge = (*clust).charge();
457  clustChargeLayer1Modules[i]->Fill(charge);
458  recHitYResLayer1Modules[i]->Fill(res_y);
459  recHitYPullLayer1Modules[i]->Fill(pull_y);
460  } else if (tTopo->pxbLayer(detId) == 2) {
461  float charge = (*clust).charge();
462  clustChargeLayer2Modules[i]->Fill(charge);
463  recHitYResLayer2Modules[i]->Fill(res_y);
464  recHitYPullLayer2Modules[i]->Fill(pull_y);
465  } else if (tTopo->pxbLayer(detId) == 3) {
466  float charge = (*clust).charge();
467  clustChargeLayer3Modules[i]->Fill(charge);
468  recHitYResLayer3Modules[i]->Fill(res_y);
469  recHitYPullLayer3Modules[i]->Fill(pull_y);
470  }
471  }
472  }
473  int sizeX = (*clust).sizeX();
474  if (tTopo->pxbLayer(detId) == 1)
475  clustXSizeLayer[0]->Fill(sizeX);
476  if (tTopo->pxbLayer(detId) == 2)
477  clustXSizeLayer[1]->Fill(sizeX);
478  if (tTopo->pxbLayer(detId) == 3)
479  clustXSizeLayer[2]->Fill(sizeX);
480 }
481 
483  const PSimHit& simHit,
484  DetId detId,
485  const PixelGeomDetUnit* theGeomDet,
486  const TrackerTopology* tTopo) {
487  int rows = theGeomDet->specificTopology().nrows();
488  int cols = theGeomDet->specificTopology().ncolumns();
489 
490  const float cmtomicron = 10000.0;
491 
492  int bunch = simHit.eventId().bunchCrossing();
493  int event = simHit.eventId().event();
494 
495  recHitBunchF->Fill(bunch);
496  if (bunch == 0)
498 
499  LocalPoint lp = recHit.localPosition();
500  float lp_x = lp.x();
501  float lp_y = lp.y();
502 
503  LocalError lerr = recHit.localPositionError();
504  float lerr_x = sqrt(lerr.xx());
505  float lerr_y = sqrt(lerr.yy());
506 
507  float sim_x1 = simHit.entryPoint().x();
508  float sim_x2 = simHit.exitPoint().x();
509  float sim_xpos = 0.5 * (sim_x1 + sim_x2);
510 
511  float sim_y1 = simHit.entryPoint().y();
512  float sim_y2 = simHit.exitPoint().y();
513  float sim_ypos = 0.5 * (sim_y1 + sim_y2);
514 
515  float pull_x = (lp_x - sim_xpos) / lerr_x;
516  float pull_y = (lp_y - sim_ypos) / lerr_y;
517 
518  if (rows == 80) {
520  } else if (rows == 160) {
522  }
523 
524  if (cols == 104) {
526  } else if (cols == 156) {
528  } else if (cols == 208) {
530  } else if (cols == 260) {
532  }
533 
534  float res_x = (lp.x() - sim_xpos) * cmtomicron;
535 
536  recHitXResAllF->Fill(res_x);
537  recHitXPullAllF->Fill(pull_x);
538 
539  float res_y = (lp.y() - sim_ypos) * cmtomicron;
540 
541  recHitYPullAllF->Fill(pull_y);
542 
543  // get cluster
544  SiPixelRecHit::ClusterRef const& clust = recHit.cluster();
545 
546  // fill plaquette dependent info
547  for (unsigned int i = 0; i < 7; i++) {
548  if (tTopo->pxfModule(detId) == i + 1) {
549  if (tTopo->pxfDisk(detId) == 1) {
550  int sizeX = (*clust).sizeX();
552 
553  int sizeY = (*clust).sizeY();
555 
556  float charge = (*clust).charge();
558 
561 
564  } else {
565  int sizeX = (*clust).sizeX();
567 
568  int sizeY = (*clust).sizeY();
570 
571  float charge = (*clust).charge();
573 
576 
579 
580  } // end else
581  } // end if module
582  else if (tTopo->pxfPanel(detId) == 2 && (tTopo->pxfModule(detId) + 4) == i + 1) {
583  if (tTopo->pxfDisk(detId) == 1) {
584  int sizeX = (*clust).sizeX();
586 
587  int sizeY = (*clust).sizeY();
589 
590  float charge = (*clust).charge();
592 
595 
598  } else {
599  int sizeX = (*clust).sizeX();
601 
602  int sizeY = (*clust).sizeY();
604 
605  float charge = (*clust).charge();
607 
610 
613 
614  } // end else
615  } // end else
616  } // end for
617 }
RunNumber_t run() const
Definition: EventID.h:38
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:81
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
MonitorElement * clustYSizeModule[8]
MonitorElement * clustXSizeDisk1Plaquettes[7]
float xx() const
Definition: LocalError.h:22
virtual int nrows() const =0
MonitorElement * recHitNsimHitDisk1
MonitorElement * recHitXPullDisk2Plaquettes[7]
MonitorElement * recHitXResFlippedLadderLayers[3]
const_iterator end(bool update=false) const
int event() const
get the contents of the subdetector field (should be protected?)
MonitorElement * recHitXFullModules
MonitorElement * recHitXResAllF
T perp() const
Definition: PV3DBase.h:69
LocalError localPositionError() const final
MonitorElement * recHitXHalfModules
MonitorElement * recHitYPullAllF
const float cmtomicron
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
MonitorElement * recHitXResDisk1Plaquettes[7]
MonitorElement * recHitEventF
MonitorElement * clustXSizeDisk2Plaquettes[7]
MonitorElement * recHitYResAllB
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
unsigned int pxfDisk(const DetId &id) const
MonitorElement * recHitYPullDisk2Plaquettes[7]
MonitorElement * recHitNsimHitLayer[3]
MonitorElement * clustXSizeLayer[3]
T y() const
Definition: PV3DBase.h:60
MonitorElement * clustChargeLayer1Modules[8]
MonitorElement * recHitYPlaquetteSize5
unsigned int pxbModule(const DetId &id) const
MonitorElement * recHitEventB
void analyze(const edm::Event &e, const edm::EventSetup &c) override
MonitorElement * recHitXPlaquetteSize1
MonitorElement * recHitXResDisk2Plaquettes[7]
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
void bookHistograms(DQMStore::IBooker &ibooker, const edm::Run &run, const edm::EventSetup &es) override
MonitorElement * recHitXPullNonFlippedLadderLayers[3]
MonitorElement * recHitYPullDisk1Plaquettes[7]
MonitorElement * recHitXPullAllB
MonitorElement * recHitXResAllB
void Fill(long long x)
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
MonitorElement * recHitBunchB
MonitorElement * clustChargeLayer2Modules[8]
TrackerHitAssociator::Config trackerHitAssociatorConfig_
MonitorElement * recHitYResDisk2Plaquettes[7]
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * recHitXResNonFlippedLadderLayers[3]
MonitorElement * recHitYResDisk1Plaquettes[7]
MonitorElement * recHitYResLayer2Modules[8]
MonitorElement * recHitYPlaquetteSize2
float yy() const
Definition: LocalError.h:24
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:19
int bunchCrossing() const
get the detector field from this detid
MonitorElement * recHitYPullLayer3Modules[8]
edm::EDGetTokenT< SiPixelRecHitCollection > siPixelRecHitCollectionToken_
MonitorElement * clustYSizeDisk1Plaquettes[7]
void beginJob() override
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
MonitorElement * recHitBunchF
void fillForward(const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
EncodedEventId eventId() const
Definition: PSimHit.h:108
unsigned int pxfModule(const DetId &id) const
MonitorElement * recHitXPlaquetteSize2
unsigned int pxbLayer(const DetId &id) const
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
MonitorElement * clustChargeLayer3Modules[8]
Definition: DetId.h:17
MonitorElement * recHitYResAllF
MonitorElement * recHitYResLayer3Modules[8]
ClusterRef cluster() const
Definition: SiPixelRecHit.h:47
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
MonitorElement * recHitYPullLayer2Modules[8]
MonitorElement * clustYSizeDisk2Plaquettes[7]
const_iterator find(id_type i, bool update=false) const
edm::EventID id() const
Definition: EventBase.h:59
MonitorElement * recHitYPlaquetteSize4
HLT enums.
virtual int ncolumns() const =0
iterator end()
Definition: DetSetNew.h:56
MonitorElement * recHitXPullFlippedLadderLayers[3]
MonitorElement * recHitNsimHitDisk2
T get() const
Definition: EventSetup.h:73
const TrackerGeomDet * idToDet(DetId) const override
MonitorElement * recHitYResLayer1Modules[8]
MonitorElement * recHitYPullLayer1Modules[8]
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
LocalPoint localPosition() const final
MonitorElement * recHitXPullDisk1Plaquettes[7]
SiPixelRecHitsValid(const edm::ParameterSet &conf)
MonitorElement * recHitYPullAllB
T x() const
Definition: PV3DBase.h:59
T const * product() const
Definition: ESHandle.h:86
MonitorElement * recHitXPullAllF
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
MonitorElement * recHitYPlaquetteSize3
MonitorElement * clustChargeDisk1Plaquettes[7]
unsigned int pxfPanel(const DetId &id) const
MonitorElement * clustChargeDisk2Plaquettes[7]
void fillBarrel(const SiPixelRecHit &, const PSimHit &, DetId, const PixelGeomDetUnit *, const TrackerTopology *tTopo)
Definition: event.py:1
Definition: Run.h:45
Our base class.
Definition: SiPixelRecHit.h:23
MonitorElement * recHitYAllModules
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
iterator begin()
Definition: DetSetNew.h:54