CMS 3D CMS Logo

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