CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelErrorEstimation.cc
Go to the documentation of this file.
1 
2 // Package: SiPixelErrorEstimation
3 // Class: SiPixelErrorEstimation
4 // Original Author: Gavril Giurgiu (JHU)
5 // Created: Fri May 4 17:48:24 CDT 2007
6 
7 #include <iostream>
9 
12 
16 
22 
24 
26 
27 using namespace std;
28 using namespace edm;
29 
30 SiPixelErrorEstimation::SiPixelErrorEstimation(const edm::ParameterSet& ps):tfile_(0), ttree_all_hits_(0), ttree_track_hits_(0)
31 {
32  //Read config file
33  outputFile_ = ps.getUntrackedParameter<string>( "outputFile", "SiPixelErrorEstimation_Ntuple.root" );
34 
35  // Replace "ctfWithMaterialTracks" with "generalTracks"
36  //src_ = ps.getUntrackedParameter<std::string>( "src", "ctfWithMaterialTracks" );
37  src_ = ps.getUntrackedParameter<std::string>( "src", "generalTracks" );
38 
39  checkType_ = ps.getParameter<bool>( "checkType" );
40  genType_ = ps.getParameter<int>( "genType" );
41  include_trk_hits_ = ps.getParameter<bool>( "include_trk_hits" );
42 }
43 
45 {}
46 
48 {
49  int bufsize = 64000;
50 
51  if ( include_trk_hits_ )
52  {
53  //tfile_ = new TFile ("SiPixelErrorEstimation_Ntuple.root" , "RECREATE");
54  //const char* tmp_name = outputFile_.c_str();
55  tfile_ = new TFile ( outputFile_.c_str() , "RECREATE");
56 
57  ttree_track_hits_ = new TTree("TrackHitNtuple", "TrackHitNtuple");
58 
59  ttree_track_hits_->Branch("evt", &evt, "evt/I", bufsize);
60  ttree_track_hits_->Branch("run", &run, "run/I", bufsize);
61 
62  ttree_track_hits_->Branch("subdetId", &subdetId, "subdetId/I", bufsize);
63 
64  ttree_track_hits_->Branch("layer" , &layer , "layer/I" , bufsize);
65  ttree_track_hits_->Branch("ladder", &ladder, "ladder/I", bufsize);
66  ttree_track_hits_->Branch("mod" , &mod , "mod/I" , bufsize);
67 
68  ttree_track_hits_->Branch("side" , &side , "side/I" , bufsize);
69  ttree_track_hits_->Branch("disk" , &disk , "disk/I" , bufsize);
70  ttree_track_hits_->Branch("blade" , &blade , "blade/I" , bufsize);
71  ttree_track_hits_->Branch("panel" , &panel , "panel/I" , bufsize);
72  ttree_track_hits_->Branch("plaq" , &plaq , "plaq/I" , bufsize);
73 
74  ttree_track_hits_->Branch("half" , &half , "half/I" , bufsize);
75  ttree_track_hits_->Branch("flipped", &flipped, "flipped/I", bufsize);
76 
77  ttree_track_hits_->Branch("rechitx", &rechitx, "rechitx/F" , bufsize);
78  ttree_track_hits_->Branch("rechity", &rechity, "rechity/F" , bufsize);
79  ttree_track_hits_->Branch("rechitz", &rechitz, "rechitz/F" , bufsize);
80 
81  ttree_track_hits_->Branch("rechiterrx", &rechiterrx, "rechiterrx/F" , bufsize);
82  ttree_track_hits_->Branch("rechiterry", &rechiterry, "rechiterry/F" , bufsize);
83 
84  ttree_track_hits_->Branch("rechitresx", &rechitresx, "rechitresx/F" , bufsize);
85  ttree_track_hits_->Branch("rechitresy", &rechitresy, "rechitresy/F" , bufsize);
86 
87  ttree_track_hits_->Branch("rechitpullx", &rechitpullx, "rechitpullx/F", bufsize);
88  ttree_track_hits_->Branch("rechitpully", &rechitpully, "rechitpully/F", bufsize);
89 
90  ttree_track_hits_->Branch("npix" , &npix , "npix/I" , bufsize);
91  ttree_track_hits_->Branch("nxpix" , &nxpix , "nxpix/I" , bufsize);
92  ttree_track_hits_->Branch("nypix" , &nypix , "nypix/I" , bufsize);
93  ttree_track_hits_->Branch("charge", &charge, "charge/F", bufsize);
94 
95  ttree_track_hits_->Branch("edgex", &edgex, "edgex/I", bufsize);
96  ttree_track_hits_->Branch("edgey", &edgey, "edgey/I", bufsize);
97 
98  ttree_track_hits_->Branch("bigx", &bigx, "bigx/I", bufsize);
99  ttree_track_hits_->Branch("bigy", &bigy, "bigy/I", bufsize);
100 
101  ttree_track_hits_->Branch("alpha", &alpha, "alpha/F", bufsize);
102  ttree_track_hits_->Branch("beta" , &beta , "beta/F" , bufsize);
103 
104  ttree_track_hits_->Branch("trk_alpha", &trk_alpha, "trk_alpha/F", bufsize);
105  ttree_track_hits_->Branch("trk_beta" , &trk_beta , "trk_beta/F" , bufsize);
106 
107  ttree_track_hits_->Branch("phi", &phi, "phi/F", bufsize);
108  ttree_track_hits_->Branch("eta", &eta, "eta/F", bufsize);
109 
110  ttree_track_hits_->Branch("simhitx", &simhitx, "simhitx/F", bufsize);
111  ttree_track_hits_->Branch("simhity", &simhity, "simhity/F", bufsize);
112 
113  ttree_track_hits_->Branch("nsimhit", &nsimhit, "nsimhit/I", bufsize);
114  ttree_track_hits_->Branch("pidhit" , &pidhit , "pidhit/I" , bufsize);
115  ttree_track_hits_->Branch("simproc", &simproc, "simproc/I", bufsize);
116 
117  ttree_track_hits_->Branch("hit_probx" , &hit_probx , "hit_probx/F" , bufsize);
118  ttree_track_hits_->Branch("hit_proby" , &hit_proby , "hit_proby/F" , bufsize);
119  ttree_track_hits_->Branch("hit_cprob0", &hit_cprob0, "hit_cprob0/F", bufsize);
120  ttree_track_hits_->Branch("hit_cprob1", &hit_cprob1, "hit_cprob1/F", bufsize);
121  ttree_track_hits_->Branch("hit_cprob2", &hit_cprob2, "hit_cprob2/F", bufsize);
122 
123  } // if ( include_trk_hits_ )
124 
125  // ----------------------------------------------------------------------
126 
127  ttree_all_hits_ = new TTree("AllHitNtuple", "AllHitNtuple");
128 
129  ttree_all_hits_->Branch("evt", &evt, "evt/I", bufsize);
130  ttree_all_hits_->Branch("run", &run, "run/I", bufsize);
131 
132  ttree_all_hits_->Branch("subdetid", &all_subdetid, "subdetid/I", bufsize);
133 
134  ttree_all_hits_->Branch("layer" , &all_layer , "layer/I" , bufsize);
135  ttree_all_hits_->Branch("ladder", &all_ladder, "ladder/I", bufsize);
136  ttree_all_hits_->Branch("mod" , &all_mod , "mod/I" , bufsize);
137 
138  ttree_all_hits_->Branch("side" , &all_side , "side/I" , bufsize);
139  ttree_all_hits_->Branch("disk" , &all_disk , "disk/I" , bufsize);
140  ttree_all_hits_->Branch("blade" , &all_blade , "blade/I" , bufsize);
141  ttree_all_hits_->Branch("panel" , &all_panel , "panel/I" , bufsize);
142  ttree_all_hits_->Branch("plaq" , &all_plaq , "plaq/I" , bufsize);
143 
144  ttree_all_hits_->Branch("half" , &all_half , "half/I" , bufsize);
145  ttree_all_hits_->Branch("flipped", &all_flipped, "flipped/I", bufsize);
146 
147  ttree_all_hits_->Branch("cols", &all_cols, "cols/I", bufsize);
148  ttree_all_hits_->Branch("rows", &all_rows, "rows/I", bufsize);
149 
150  ttree_all_hits_->Branch("rechitx" , &all_rechitx , "rechitx/F" , bufsize);
151  ttree_all_hits_->Branch("rechity" , &all_rechity , "rechity/F" , bufsize);
152  ttree_all_hits_->Branch("rechitz" , &all_rechitz , "rechitz/F" , bufsize);
153 
154  ttree_all_hits_->Branch("rechiterrx" , &all_rechiterrx , "rechiterrx/F" , bufsize);
155  ttree_all_hits_->Branch("rechiterry" , &all_rechiterry , "rechiterry/F" , bufsize);
156 
157  ttree_all_hits_->Branch("rechitresx" , &all_rechitresx , "rechitresx/F" , bufsize);
158  ttree_all_hits_->Branch("rechitresy" , &all_rechitresy , "rechitresy/F" , bufsize);
159 
160  ttree_all_hits_->Branch("rechitpullx", &all_rechitpullx, "rechitpullx/F", bufsize);
161  ttree_all_hits_->Branch("rechitpully", &all_rechitpully, "rechitpully/F", bufsize);
162 
163  ttree_all_hits_->Branch("npix" , &all_npix , "npix/I" , bufsize);
164  ttree_all_hits_->Branch("nxpix" , &all_nxpix , "nxpix/I" , bufsize);
165  ttree_all_hits_->Branch("nypix" , &all_nypix , "nypix/I" , bufsize);
166 
167  ttree_all_hits_->Branch("edgex", &all_edgex, "edgex/I", bufsize);
168  ttree_all_hits_->Branch("edgey", &all_edgey, "edgey/I", bufsize);
169 
170  ttree_all_hits_->Branch("bigx", &all_bigx, "bigx/I", bufsize);
171  ttree_all_hits_->Branch("bigy", &all_bigy, "bigy/I", bufsize);
172 
173  ttree_all_hits_->Branch("alpha", &all_alpha, "alpha/F", bufsize);
174  ttree_all_hits_->Branch("beta" , &all_beta , "beta/F" , bufsize);
175 
176  ttree_all_hits_->Branch("simhitx", &all_simhitx, "simhitx/F", bufsize);
177  ttree_all_hits_->Branch("simhity", &all_simhity, "simhity/F", bufsize);
178 
179  ttree_all_hits_->Branch("nsimhit", &all_nsimhit, "nsimhit/I", bufsize);
180  ttree_all_hits_->Branch("pidhit" , &all_pidhit , "pidhit/I" , bufsize);
181  ttree_all_hits_->Branch("simproc", &all_simproc, "simproc/I", bufsize);
182 
183  ttree_all_hits_->Branch("vtxr", &all_vtxr, "vtxr/F", bufsize);
184  ttree_all_hits_->Branch("vtxz", &all_vtxz, "vtxz/F", bufsize);
185 
186  ttree_all_hits_->Branch("simpx", &all_simpx, "simpx/F", bufsize);
187  ttree_all_hits_->Branch("simpy", &all_simpy, "simpy/F", bufsize);
188  ttree_all_hits_->Branch("simpz", &all_simpz, "simpz/F", bufsize);
189 
190  ttree_all_hits_->Branch("eloss", &all_eloss, "eloss/F", bufsize);
191 
192  ttree_all_hits_->Branch("simphi", &all_simphi, "simphi/F", bufsize);
193  ttree_all_hits_->Branch("simtheta", &all_simtheta, "simtheta/F", bufsize);
194 
195  ttree_all_hits_->Branch("trkid", &all_trkid, "trkid/I", bufsize);
196 
197  ttree_all_hits_->Branch("x1", &all_x1, "x1/F", bufsize);
198  ttree_all_hits_->Branch("x2", &all_x2, "x2/F", bufsize);
199  ttree_all_hits_->Branch("y1", &all_x1, "y1/F", bufsize);
200  ttree_all_hits_->Branch("y2", &all_x2, "y2/F", bufsize);
201  ttree_all_hits_->Branch("z1", &all_x1, "z1/F", bufsize);
202  ttree_all_hits_->Branch("z2", &all_x2, "z2/F", bufsize);
203 
204  ttree_all_hits_->Branch("row1", &all_row1, "row1/F", bufsize);
205  ttree_all_hits_->Branch("row2", &all_row2, "row2/F", bufsize);
206  ttree_all_hits_->Branch("col1", &all_col1, "col1/F", bufsize);
207  ttree_all_hits_->Branch("col2", &all_col2, "col2/F", bufsize);
208 
209  ttree_all_hits_->Branch("gx1", &all_gx1, "gx1/F", bufsize);
210  ttree_all_hits_->Branch("gx2", &all_gx2, "gx2/F", bufsize);
211  ttree_all_hits_->Branch("gy1", &all_gx1, "gy1/F", bufsize);
212  ttree_all_hits_->Branch("gy2", &all_gx2, "gy2/F", bufsize);
213  ttree_all_hits_->Branch("gz1", &all_gx1, "gz1/F", bufsize);
214  ttree_all_hits_->Branch("gz2", &all_gx2, "gz2/F", bufsize);
215 
216  ttree_all_hits_->Branch("simtrketa", &all_simtrketa, "simtrketa/F", bufsize);
217  ttree_all_hits_->Branch("simtrkphi", &all_simtrkphi, "simtrkphi/F", bufsize);
218 
219  ttree_all_hits_->Branch("clust_row", &all_clust_row, "clust_row/F", bufsize);
220  ttree_all_hits_->Branch("clust_col", &all_clust_col, "clust_col/F", bufsize);
221 
222  ttree_all_hits_->Branch("clust_x", &all_clust_x, "clust_x/F", bufsize);
223  ttree_all_hits_->Branch("clust_y", &all_clust_y, "clust_y/F", bufsize);
224 
225  ttree_all_hits_->Branch("clust_q", &all_clust_q, "clust_q/F", bufsize);
226 
227  ttree_all_hits_->Branch("clust_maxpixcol", &all_clust_maxpixcol, "clust_maxpixcol/I", bufsize);
228  ttree_all_hits_->Branch("clust_maxpixrow", &all_clust_maxpixrow, "clust_maxpixrow/I", bufsize);
229  ttree_all_hits_->Branch("clust_minpixcol", &all_clust_minpixcol, "clust_minpixcol/I", bufsize);
230  ttree_all_hits_->Branch("clust_minpixrow", &all_clust_minpixrow, "clust_minpixrow/I", bufsize);
231 
232  ttree_all_hits_->Branch("clust_geoid", &all_clust_geoid, "clust_geoid/I", bufsize);
233 
234  ttree_all_hits_->Branch("clust_alpha", &all_clust_alpha, "clust_alpha/F", bufsize);
235  ttree_all_hits_->Branch("clust_beta" , &all_clust_beta , "clust_beta/F" , bufsize);
236 
237  ttree_all_hits_->Branch("rowpix", all_pixrow, "row[npix]/F", bufsize);
238  ttree_all_hits_->Branch("colpix", all_pixcol, "col[npix]/F", bufsize);
239  ttree_all_hits_->Branch("adc", all_pixadc, "adc[npix]/F", bufsize);
240  ttree_all_hits_->Branch("xpix", all_pixx, "x[npix]/F", bufsize);
241  ttree_all_hits_->Branch("ypix", all_pixy, "y[npix]/F", bufsize);
242  ttree_all_hits_->Branch("gxpix", all_pixgx, "gx[npix]/F", bufsize);
243  ttree_all_hits_->Branch("gypix", all_pixgy, "gy[npix]/F", bufsize);
244  ttree_all_hits_->Branch("gzpix", all_pixgz, "gz[npix]/F", bufsize);
245 
246  ttree_all_hits_->Branch("all_hit_probx", &all_hit_probx, "all_hit_probx/F" , bufsize);
247  ttree_all_hits_->Branch("all_hit_proby", &all_hit_proby, "all_hit_proby/F" , bufsize);
248  ttree_all_hits_->Branch("all_hit_cprob0", &all_hit_cprob0, "all_hit_cprob0/F", bufsize);
249  ttree_all_hits_->Branch("all_hit_cprob1", &all_hit_cprob1, "all_hit_cprob1/F", bufsize);
250  ttree_all_hits_->Branch("all_hit_cprob2", &all_hit_cprob2, "all_hit_cprob2/F", bufsize);
251 
252 }
253 
255 {
256  tfile_->Write();
257  tfile_->Close();
258 }
259 
260 void
262 {
263  using namespace edm;
264 
265  run = e.id().run();
266  evt = e.id().event();
267 
268  if ( evt%1000 == 0 )
269  cout << "evt = " << evt << endl;
270 
271  float math_pi = 3.14159265;
272  float radtodeg = 180.0 / math_pi;
273 
274  DetId detId;
275 
278  float mindist = 999999.9;
279 
280  std::vector<PSimHit> matched;
281  TrackerHitAssociator associate(e);
282 
284  es.get<TrackerDigiGeometryRecord> ().get (pDD);
285  const TrackerGeometry* tracker = &(* pDD);
286 
287  // --------------------------------------- all hits -----------------------------------------------------------
289  e.getByLabel( "siPixelRecHits", recHitColl);
290 
292  e.getByLabel("g4SimHits", simtracks);
293 
294  //-----Iterate over detunits
295  for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++)
296  {
297  DetId detId = ((*it)->geographicalId());
298 
299  SiPixelRecHitCollection::const_iterator dsmatch = recHitColl->find(detId);
300  if (dsmatch == recHitColl->end()) continue;
301 
302  SiPixelRecHitCollection::DetSet pixelrechitRange = *dsmatch;
303  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorBegin = pixelrechitRange.begin();
304  SiPixelRecHitCollection::DetSet::const_iterator pixelrechitRangeIteratorEnd = pixelrechitRange.end();
305  SiPixelRecHitCollection::DetSet::const_iterator pixeliter = pixelrechitRangeIteratorBegin;
306  std::vector<PSimHit> matched;
307 
308  //----Loop over rechits for this detId
309  for ( ; pixeliter != pixelrechitRangeIteratorEnd; ++pixeliter)
310  {
311  matched.clear();
312  matched = associate.associateHit(*pixeliter);
313  // only consider rechits that have associated simhit
314  // otherwise cannot determine residiual
315  if ( matched.empty() )
316  {
317  cout << "SiPixelErrorEstimation::analyze: rechits without associated simhit !!!!!!!"
318  << endl;
319  continue;
320  }
321 
322  all_subdetid = -9999;
323 
324  all_layer = -9999;
325  all_ladder = -9999;
326  all_mod = -9999;
327 
328  all_side = -9999;
329  all_disk = -9999;
330  all_blade = -9999;
331  all_panel = -9999;
332  all_plaq = -9999;
333 
334  all_half = -9999;
335  all_flipped = -9999;
336 
337  all_cols = -9999;
338  all_rows = -9999;
339 
340  all_rechitx = -9999;
341  all_rechity = -9999;
342  all_rechitz = -9999;
343 
344  all_simhitx = -9999;
345  all_simhity = -9999;
346 
347  all_rechiterrx = -9999;
348  all_rechiterry = -9999;
349 
350  all_rechitresx = -9999;
351  all_rechitresy = -9999;
352 
353  all_rechitpullx = -9999;
354  all_rechitpully = -9999;
355 
356  all_npix = -9999;
357  all_nxpix = -9999;
358  all_nypix = -9999;
359 
360  all_edgex = -9999;
361  all_edgey = -9999;
362 
363  all_bigx = -9999;
364  all_bigy = -9999;
365 
366  all_alpha = -9999;
367  all_beta = -9999;
368 
369  all_simphi = -9999;
370  all_simtheta = -9999;
371 
372  all_simhitx = -9999;
373  all_simhity = -9999;
374 
375  all_nsimhit = -9999;
376  all_pidhit = -9999;
377  all_simproc = -9999;
378 
379  all_vtxr = -9999;
380  all_vtxz = -9999;
381 
382  all_simpx = -9999;
383  all_simpy = -9999;
384  all_simpz = -9999;
385 
386  all_eloss = -9999;
387 
388  all_trkid = -9999;
389 
390  all_x1 = -9999;
391  all_x2 = -9999;
392  all_y1 = -9999;
393  all_y2 = -9999;
394  all_z1 = -9999;
395  all_z2 = -9999;
396 
397  all_row1 = -9999;
398  all_row2 = -9999;
399  all_col1 = -9999;
400  all_col2 = -9999;
401 
402  all_gx1 = -9999;
403  all_gx2 = -9999;
404  all_gy1 = -9999;
405  all_gy2 = -9999;
406  all_gz1 = -9999;
407  all_gz2 = -9999;
408 
409  all_simtrketa = -9999;
410  all_simtrkphi = -9999;
411 
412  all_clust_row = -9999;
413  all_clust_col = -9999;
414 
415  all_clust_x = -9999;
416  all_clust_y = -9999;
417 
418  all_clust_q = -9999;
419 
420  all_clust_maxpixcol = -9999;
421  all_clust_maxpixrow = -9999;
422  all_clust_minpixcol = -9999;
423  all_clust_minpixrow = -9999;
424 
425  all_clust_geoid = -9999;
426 
427  all_clust_alpha = -9999;
428  all_clust_beta = -9999;
429 
430  /*
431  for (int i=0; i<all_npix; ++i)
432  {
433  all_pixrow[i] = -9999;
434  all_pixcol[i] = -9999;
435  all_pixadc[i] = -9999;
436  all_pixx[i] = -9999;
437  all_pixy[i] = -9999;
438  all_pixgx[i] = -9999;
439  all_pixgy[i] = -9999;
440  all_pixgz[i] = -9999;
441  }
442  */
443 
444  all_hit_probx = -9999;
445  all_hit_proby = -9999;
446  all_hit_cprob0 = -9999;
447  all_hit_cprob1 = -9999;
448  all_hit_cprob2 = -9999;
449 
450  all_nsimhit = (int)matched.size();
451 
452  all_subdetid = (int)detId.subdetId();
453  // only consider rechits in pixel barrel and pixel forward
454  if ( !(all_subdetid==1 || all_subdetid==2) )
455  {
456  cout << "SiPixelErrorEstimation::analyze: Not in a pixel detector !!!!!" << endl;
457  continue;
458  }
459 
460  const PixelGeomDetUnit* theGeomDet
461  = dynamic_cast<const PixelGeomDetUnit*> ( tracker->idToDet(detId) );
462 
463  const PixelTopology* topol = &(theGeomDet->specificTopology());
464 
465  const int maxPixelCol = pixeliter->cluster()->maxPixelCol();
466  const int maxPixelRow = pixeliter->cluster()->maxPixelRow();
467  const int minPixelCol = pixeliter->cluster()->minPixelCol();
468  const int minPixelRow = pixeliter->cluster()->minPixelRow();
469 
470  //all_hit_probx = (float)pixeliter->probabilityX();
471  //all_hit_proby = (float)pixeliter->probabilityY();
472  all_hit_cprob0 = (float)pixeliter->clusterProbability(0);
473  all_hit_cprob1 = (float)pixeliter->clusterProbability(1);
474  all_hit_cprob2 = (float)pixeliter->clusterProbability(2);
475 
476  // check whether the cluster is at the module edge
477  if ( topol->isItEdgePixelInX( minPixelRow ) ||
478  topol->isItEdgePixelInX( maxPixelRow ) )
479  all_edgex = 1;
480  else
481  all_edgex = 0;
482 
483  if ( topol->isItEdgePixelInY( minPixelCol ) ||
484  topol->isItEdgePixelInY( maxPixelCol ) )
485  all_edgey = 1;
486  else
487  all_edgey = 0;
488 
489  // check whether this rechit contains big pixels
490  if ( topol->containsBigPixelInX(minPixelRow, maxPixelRow) )
491  all_bigx = 1;
492  else
493  all_bigx = 0;
494 
495  if ( topol->containsBigPixelInY(minPixelCol, maxPixelCol) )
496  all_bigy = 1;
497  else
498  all_bigy = 0;
499 
500  if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
501  {
502  PXBDetId bdetid(detId);
503  all_layer = bdetid.layer();
504  all_ladder = bdetid.ladder();
505  all_mod = bdetid.module();
506 
507  int tmp_nrows = theGeomDet->specificTopology().nrows();
508  if ( tmp_nrows == 80 )
509  all_half = 1;
510  else if ( tmp_nrows == 160 )
511  all_half = 0;
512  else
513  cout << "-------------------------------------------------- Wrong module size !!!" << endl;
514 
515  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
516  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
517 
518  if ( tmp2<tmp1 )
519  all_flipped = 1;
520  else
521  all_flipped = 0;
522  }
523  else if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelEndcap )
524  {
525  PXFDetId fdetid(detId);
526  all_side = fdetid.side();
527  all_disk = fdetid.disk();
528  all_blade = fdetid.blade();
529  all_panel = fdetid.panel();
530  all_plaq = fdetid.module(); // also known as plaquette
531 
532  } // else if ( detId.subdetId()==PixelSubdetector::PixelEndcap )
533  else std::cout << "We are not in the pixel detector" << (int)detId.subdetId() << endl;
534 
535  all_cols = theGeomDet->specificTopology().ncolumns();
536  all_rows = theGeomDet->specificTopology().nrows();
537 
538  LocalPoint lp = pixeliter->localPosition();
539  // gavril: change this name
540  all_rechitx = lp.x();
541  all_rechity = lp.y();
542  all_rechitz = lp.z();
543 
544  LocalError le = pixeliter->localPositionError();
545  all_rechiterrx = sqrt( le.xx() );
546  all_rechiterry = sqrt( le.yy() );
547 
548  bool found_hit_from_generated_particle = false;
549 
550  //---Loop over sim hits, fill closest
551  float closest_dist = 99999.9;
552  std::vector<PSimHit>::const_iterator closest_simhit = matched.begin();
553 
554  for (std::vector<PSimHit>::const_iterator m = matched.begin(); m < matched.end(); m++)
555  {
556  if ( checkType_ )
557  {
558  int pid = (*m).particleType();
559  if ( abs(pid) != genType_ )
560  continue;
561  }
562 
563  float simhitx = 0.5 * ( (*m).entryPoint().x() + (*m).exitPoint().x() );
564  float simhity = 0.5 * ( (*m).entryPoint().y() + (*m).exitPoint().y() );
565 
566  float x_res = simhitx - rechitx;
567  float y_res = simhity - rechity;
568 
569  float dist = sqrt( x_res*x_res + y_res*y_res );
570 
571  if ( dist < closest_dist )
572  {
573  closest_dist = dist;
574  closest_simhit = m;
575  found_hit_from_generated_particle = true;
576  }
577  } // end sim hit loop
578 
579  // If this recHit does not have any simHit with the same particleType as the particles generated
580  // ignore it as most probably comes from delta rays.
581  if ( checkType_ && !found_hit_from_generated_particle )
582  continue;
583 
584  all_x1 = (*closest_simhit).entryPoint().x(); // width (row index, in col direction)
585  all_y1 = (*closest_simhit).entryPoint().y(); // length (col index, in row direction)
586  all_z1 = (*closest_simhit).entryPoint().z();
587  all_x2 = (*closest_simhit).exitPoint().x();
588  all_y2 = (*closest_simhit).exitPoint().y();
589  all_z2 = (*closest_simhit).exitPoint().z();
590  GlobalPoint GP1 =
591  theGeomDet->surface().toGlobal( Local3DPoint( (*closest_simhit).entryPoint().x(),
592  (*closest_simhit).entryPoint().y(),
593  (*closest_simhit).entryPoint().z() ) );
594  GlobalPoint GP2 =
595  theGeomDet->surface().toGlobal (Local3DPoint( (*closest_simhit).exitPoint().x(),
596  (*closest_simhit).exitPoint().y(),
597  (*closest_simhit).exitPoint().z() ) );
598  all_gx1 = GP1.x();
599  all_gx2 = GP2.x();
600  all_gy1 = GP1.y();
601  all_gy2 = GP2.y();
602  all_gz1 = GP1.z();
603  all_gz2 = GP2.z();
604 
605  MeasurementPoint mp1 =
606  topol->measurementPosition( LocalPoint( (*closest_simhit).entryPoint().x(),
607  (*closest_simhit).entryPoint().y(),
608  (*closest_simhit).entryPoint().z() ) );
609  MeasurementPoint mp2 =
610  topol->measurementPosition( LocalPoint( (*closest_simhit).exitPoint().x(),
611  (*closest_simhit).exitPoint().y(),
612  (*closest_simhit).exitPoint().z() ) );
613  all_row1 = mp1.x();
614  all_col1 = mp1.y();
615  all_row2 = mp2.x();
616  all_col2 = mp2.y();
617 
618  all_simhitx = 0.5*(all_x1+all_x2);
619  all_simhity = 0.5*(all_y1+all_y2);
620 
623 
626 
627  SiPixelRecHit::ClusterRef const& clust = pixeliter->cluster();
628 
629  all_npix = clust->size();
630  all_nxpix = clust->sizeX();
631  all_nypix = clust->sizeY();
632 
633  all_clust_row = clust->x();
634  all_clust_col = clust->y();
635 
637  all_clust_x = lp2.x();
638  all_clust_y = lp2.y();
639 
640  all_clust_q = clust->charge();
641 
642  all_clust_maxpixcol = clust->maxPixelCol();
643  all_clust_maxpixrow = clust->maxPixelRow();
644  all_clust_minpixcol = clust->minPixelCol();
645  all_clust_minpixrow = clust->minPixelRow();
646 
647  all_clust_geoid = clust->geographicalId();
648 
649  all_simpx = (*closest_simhit).momentumAtEntry().x();
650  all_simpy = (*closest_simhit).momentumAtEntry().y();
651  all_simpz = (*closest_simhit).momentumAtEntry().z();
652  all_eloss = (*closest_simhit).energyLoss();
653  all_simphi = (*closest_simhit).phiAtEntry();
654  all_simtheta = (*closest_simhit).thetaAtEntry();
655  all_pidhit = (*closest_simhit).particleType();
656  all_trkid = (*closest_simhit).trackId();
657 
658  //--- Fill alpha and beta -- more useful for exploring the residuals...
659  all_beta = atan2(all_simpz, all_simpy);
660  all_alpha = atan2(all_simpz, all_simpx);
661 
662  all_simproc = (int)closest_simhit->processType();
663 
664  const edm::SimTrackContainer& trks = *(simtracks.product());
665  SimTrackContainer::const_iterator trksiter;
666  for (trksiter = trks.begin(); trksiter != trks.end(); trksiter++)
667  if ( (int)trksiter->trackId() == all_trkid )
668  {
669  all_simtrketa = trksiter->momentum().eta();
670  all_simtrkphi = trksiter->momentum().phi();
671  }
672 
673  all_vtxz = theGeomDet->surface().position().z();
674  all_vtxr = theGeomDet->surface().position().perp();
675 
676  //computeAnglesFromDetPosition(clust,
677  // theGeomDet,
678  // all_clust_alpha, all_clust_beta )
679 
680  const std::vector<SiPixelCluster::Pixel>& pixvector = clust->pixels();
681  for ( int i=0; i<(int)pixvector.size(); ++i)
682  {
683  SiPixelCluster::Pixel holdpix = pixvector[i];
684  all_pixrow[i] = holdpix.x;
685  all_pixcol[i] = holdpix.y;
686  all_pixadc[i] = holdpix.adc;
687  LocalPoint lp = topol->localPosition(MeasurementPoint(holdpix.x, holdpix.y));
688  all_pixx[i] = lp.x();
689  all_pixy[i]= lp.y();
690  GlobalPoint GP = theGeomDet->surface().toGlobal(Local3DPoint(lp.x(),lp.y(),lp.z()));
691  all_pixgx[i] = GP.x();
692  all_pixgy[i]= GP.y();
693  all_pixgz[i]= GP.z();
694  }
695 
696  ttree_all_hits_->Fill();
697 
698  } // for ( ; pixeliter != pixelrechitRangeIteratorEnd; ++pixeliter)
699 
700  } // for (TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++)
701 
702  // ------------------------------------------------ all hits ---------------------------------------------------------------
703 
704 
705 
706  // ------------------------------------------------ track hits only --------------------------------------------------------
707 
708  if ( include_trk_hits_ )
709  {
710  // Get tracks
711  edm::Handle<reco::TrackCollection> trackCollection;
712  e.getByLabel(src_, trackCollection);
713  const reco::TrackCollection *tracks = trackCollection.product();
714  reco::TrackCollection::const_iterator tciter;
715 
716  if ( tracks->size() > 0 )
717  {
718  // Loop on tracks
719  for ( tciter=tracks->begin(); tciter!=tracks->end(); ++tciter)
720  {
721  // First loop on hits: find matched hits
722  for ( trackingRecHit_iterator it = tciter->recHitsBegin(); it != tciter->recHitsEnd(); ++it)
723  {
724  const TrackingRecHit &thit = **it;
725  // Is it a matched hit?
726  const SiPixelRecHit* matchedhit = dynamic_cast<const SiPixelRecHit*>(&thit);
727 
728  if ( matchedhit )
729  {
730  rechitx = -9999.9;
731  rechity = -9999.9;
732  rechitz = -9999.9;
733  rechiterrx = -9999.9;
734  rechiterry = -9999.9;
735  rechitresx = -9999.9;
736  rechitresy = -9999.9;
737  rechitpullx = -9999.9;
738  rechitpully = -9999.9;
739 
740  npix = -9999;
741  nxpix = -9999;
742  nypix = -9999;
743  charge = -9999.9;
744 
745  edgex = -9999;
746  edgey = -9999;
747 
748  bigx = -9999;
749  bigy = -9999;
750 
751  alpha = -9999.9;
752  beta = -9999.9;
753 
754  phi = -9999.9;
755  eta = -9999.9;
756 
757  subdetId = -9999;
758 
759  layer = -9999;
760  ladder = -9999;
761  mod = -9999;
762  side = -9999;
763  disk = -9999;
764  blade = -9999;
765  panel = -9999;
766  plaq = -9999;
767 
768  half = -9999;
769  flipped = -9999;
770 
771  nsimhit = -9999;
772  pidhit = -9999;
773  simproc = -9999;
774 
775  simhitx = -9999.9;
776  simhity = -9999.9;
777 
778  hit_probx = -9999.9;
779  hit_proby = -9999.9;
780  hit_cprob0 = -9999.9;
781  hit_cprob1 = -9999.9;
782  hit_cprob2 = -9999.9;
783 
784  position = (*it)->localPosition();
785  error = (*it)->localPositionError();
786 
787  rechitx = position.x();
788  rechity = position.y();
789  rechitz = position.z();
790  rechiterrx = sqrt(error.xx());
791  rechiterry = sqrt(error.yy());
792 
793  npix = matchedhit->cluster()->size();
794  nxpix = matchedhit->cluster()->sizeX();
795  nypix = matchedhit->cluster()->sizeY();
796  charge = matchedhit->cluster()->charge();
797 
798  //hit_probx = (float)matchedhit->probabilityX();
799  //hit_proby = (float)matchedhit->probabilityY();
800  hit_cprob0 = (float)matchedhit->clusterProbability(0);
801  hit_cprob1 = (float)matchedhit->clusterProbability(1);
802  hit_cprob2 = (float)matchedhit->clusterProbability(2);
803 
804 
805  //Association of the rechit to the simhit
806  matched.clear();
807  matched = associate.associateHit(*matchedhit);
808 
809  nsimhit = (int)matched.size();
810 
811  if ( !matched.empty() )
812  {
813  mindist = 999999.9;
814  float distx, disty, dist;
815  bool found_hit_from_generated_particle = false;
816 
817  int n_assoc_muon = 0;
818 
819  vector<PSimHit>::const_iterator closestit = matched.begin();
820  for (vector<PSimHit>::const_iterator m=matched.begin(); m<matched.end(); ++m)
821  {
822  if ( checkType_ )
823  { // only consider associated simhits with the generated pid (muons)
824  int pid = (*m).particleType();
825  if ( abs(pid) != genType_ )
826  continue;
827  }
828 
829  float simhitx = 0.5 * ( (*m).entryPoint().x() + (*m).exitPoint().x() );
830  float simhity = 0.5 * ( (*m).entryPoint().y() + (*m).exitPoint().y() );
831 
832  distx = fabs(rechitx - simhitx);
833  disty = fabs(rechity - simhity);
834  dist = sqrt( distx*distx + disty*disty );
835 
836  if ( dist < mindist )
837  {
838  n_assoc_muon++;
839 
840  mindist = dist;
841  closestit = m;
842  found_hit_from_generated_particle = true;
843  }
844  } // for (vector<PSimHit>::const_iterator m=matched.begin(); m<matched.end(); m++)
845 
846  // This recHit does not have any simHit with the same particleType as the particles generated
847  // Ignore it as most probably come from delta rays.
848  if ( checkType_ && !found_hit_from_generated_particle )
849  continue;
850 
851  //if ( n_assoc_muon > 1 )
852  //{
853  // // cout << " ----- This is not good: n_assoc_muon = " << n_assoc_muon << endl;
854  // // cout << "evt = " << evt << endl;
855  //}
856 
857  detId = (*it)->geographicalId();
858 
859  const PixelGeomDetUnit* theGeomDet =
860  dynamic_cast<const PixelGeomDetUnit*> ((*tracker).idToDet(detId) );
861 
862  const PixelTopology* theTopol = &(theGeomDet->specificTopology());
863 
864  pidhit = (*closestit).particleType();
865  simproc = (int)(*closestit).processType();
866 
867  simhitx = 0.5*( (*closestit).entryPoint().x() + (*closestit).exitPoint().x() );
868  simhity = 0.5*( (*closestit).entryPoint().y() + (*closestit).exitPoint().y() );
869 
872  rechitpullx = ( rechitx - simhitx ) / sqrt(error.xx());
873  rechitpully = ( rechity - simhity ) / sqrt(error.yy());
874 
875  float simhitpx = (*closestit).momentumAtEntry().x();
876  float simhitpy = (*closestit).momentumAtEntry().y();
877  float simhitpz = (*closestit).momentumAtEntry().z();
878 
879  beta = atan2(simhitpz, simhitpy) * radtodeg;
880  alpha = atan2(simhitpz, simhitpx) * radtodeg;
881 
882  //beta = fabs(atan2(simhitpz, simhitpy)) * radtodeg;
883  //alpha = fabs(atan2(simhitpz, simhitpx)) * radtodeg;
884 
885  // calculate alpha and beta exactly as in PixelCPEBase.cc
886  float locx = simhitpx;
887  float locy = simhitpy;
888  float locz = simhitpz;
889 
890  bool isFlipped = false;
891  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
892  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
893  if ( tmp2<tmp1 )
894  isFlipped = true;
895  else
896  isFlipped = false;
897 
898  trk_alpha = acos(locx/sqrt(locx*locx+locz*locz)) * radtodeg;
899  if ( isFlipped ) // &&& check for FPIX !!!
900  trk_alpha = 180.0 - trk_alpha ;
901 
902  trk_beta = acos(locy/sqrt(locy*locy+locz*locz)) * radtodeg;
903 
904 
905  phi = tciter->momentum().phi() / math_pi*180.0;
906  eta = tciter->momentum().eta();
907 
908  const int maxPixelCol = (*matchedhit).cluster()->maxPixelCol();
909  const int maxPixelRow = (*matchedhit).cluster()->maxPixelRow();
910  const int minPixelCol = (*matchedhit).cluster()->minPixelCol();
911  const int minPixelRow = (*matchedhit).cluster()->minPixelRow();
912 
913  // check whether the cluster is at the module edge
914  if ( theTopol->isItEdgePixelInX( minPixelRow ) ||
915  theTopol->isItEdgePixelInX( maxPixelRow ) )
916  edgex = 1;
917  else
918  edgex = 0;
919 
920  if ( theTopol->isItEdgePixelInY( minPixelCol ) ||
921  theTopol->isItEdgePixelInY( maxPixelCol ) )
922  edgey = 1;
923  else
924  edgey = 0;
925 
926  // check whether this rechit contains big pixels
927  if ( theTopol->containsBigPixelInX(minPixelRow, maxPixelRow) )
928  bigx = 1;
929  else
930  bigx = 0;
931 
932  if ( theTopol->containsBigPixelInY(minPixelCol, maxPixelCol) )
933  bigy = 1;
934  else
935  bigy = 0;
936 
937  subdetId = (int)detId.subdetId();
938 
939  if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelBarrel )
940  {
941 
942  int tmp_nrows = theGeomDet->specificTopology().nrows();
943  if ( tmp_nrows == 80 )
944  half = 1;
945  else if ( tmp_nrows == 160 )
946  half = 0;
947  else
948  cout << "-------------------------------------------------- Wrong module size !!!" << endl;
949 
950  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
951  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
952 
953  if ( tmp2<tmp1 )
954  flipped = 1;
955  else
956  flipped = 0;
957 
958  PXBDetId bdetid(detId);
959  layer = bdetid.layer(); // Layer: 1,2,3.
960  ladder = bdetid.ladder(); // Ladder: 1-20, 32, 44.
961  mod = bdetid.module(); // Mod: 1-8.
962  }
963  else if ( (int)detId.subdetId() == (int)PixelSubdetector::PixelEndcap )
964  {
965  PXFDetId fdetid(detId);
966  side = fdetid.side();
967  disk = fdetid.disk();
968  blade = fdetid.blade();
969  panel = fdetid.panel();
970  plaq = fdetid.module(); // also known as plaquette
971 
972  float tmp1 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp();
973  float tmp2 = theGeomDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp();
974 
975  if ( tmp2<tmp1 )
976  flipped = 1;
977  else
978  flipped = 0;
979 
980  } // else if ( detId.subdetId()==PixelSubdetector::PixelEndcap )
981  //else std::// cout << "We are not in the pixel detector. detId.subdetId() = " << (int)detId.subdetId() << endl;
982 
983  ttree_track_hits_->Fill();
984 
985  } // if ( !matched.empty() )
986  else
987  cout << "---------------- RecHit with no associated SimHit !!! -------------------------- " << endl;
988 
989  } // if ( matchedhit )
990 
991  } // end of loop on hits
992 
993  } //end of loop on track
994 
995  } // tracks > 0.
996 
997  } // if ( include_trk_hits_ )
998 
999  // ----------------------------------------------- track hits only -----------------------------------------------------------
1000 
1001 }
1002 
1005  const GeomDetUnit & det,
1006  float& alpha, float& beta )
1007 {
1008  //--- This is a new det unit, so cache it
1009  const PixelGeomDetUnit* theDet = dynamic_cast<const PixelGeomDetUnit*>( &det );
1010  if (! theDet)
1011  {
1012  cout << "---------------------------------------------- Not a pixel detector !!!!!!!!!!!!!!" << endl;
1013  assert(0);
1014  }
1015 
1016  const PixelTopology* theTopol = &(theDet->specificTopology());
1017 
1018  // get cluster center of gravity (of charge)
1019  float xcenter = cl.x();
1020  float ycenter = cl.y();
1021 
1022  // get the cluster position in local coordinates (cm)
1023  LocalPoint lp = theTopol->localPosition( MeasurementPoint(xcenter, ycenter) );
1024  //float lp_mod = sqrt( lp.x()*lp.x() + lp.y()*lp.y() + lp.z()*lp.z() );
1025 
1026  // get the cluster position in global coordinates (cm)
1027  GlobalPoint gp = theDet->surface().toGlobal( lp );
1028  float gp_mod = sqrt( gp.x()*gp.x() + gp.y()*gp.y() + gp.z()*gp.z() );
1029 
1030  // normalize
1031  float gpx = gp.x()/gp_mod;
1032  float gpy = gp.y()/gp_mod;
1033  float gpz = gp.z()/gp_mod;
1034 
1035  // make a global vector out of the global point; this vector will point from the
1036  // origin of the detector to the cluster
1037  GlobalVector gv(gpx, gpy, gpz);
1038 
1039  // make local unit vector along local X axis
1040  const Local3DVector lvx(1.0, 0.0, 0.0);
1041 
1042  // get the unit X vector in global coordinates/
1043  GlobalVector gvx = theDet->surface().toGlobal( lvx );
1044 
1045  // make local unit vector along local Y axis
1046  const Local3DVector lvy(0.0, 1.0, 0.0);
1047 
1048  // get the unit Y vector in global coordinates
1049  GlobalVector gvy = theDet->surface().toGlobal( lvy );
1050 
1051  // make local unit vector along local Z axis
1052  const Local3DVector lvz(0.0, 0.0, 1.0);
1053 
1054  // get the unit Z vector in global coordinates
1055  GlobalVector gvz = theDet->surface().toGlobal( lvz );
1056 
1057  // calculate the components of gv (the unit vector pointing to the cluster)
1058  // in the local coordinate system given by the basis {gvx, gvy, gvz}
1059  // note that both gv and the basis {gvx, gvy, gvz} are given in global coordinates
1060  float gv_dot_gvx = gv.x()*gvx.x() + gv.y()*gvx.y() + gv.z()*gvx.z();
1061  float gv_dot_gvy = gv.x()*gvy.x() + gv.y()*gvy.y() + gv.z()*gvy.z();
1062  float gv_dot_gvz = gv.x()*gvz.x() + gv.y()*gvz.y() + gv.z()*gvz.z();
1063 
1064  // calculate angles
1065  alpha = atan2( gv_dot_gvz, gv_dot_gvx );
1066  beta = atan2( gv_dot_gvz, gv_dot_gvy );
1067 
1068  // calculate cotalpha and cotbeta
1069  // cotalpha_ = 1.0/tan(alpha_);
1070  // cotbeta_ = 1.0/tan(beta_ );
1071  // or like this
1072  //cotalpha_ = gv_dot_gvx / gv_dot_gvz;
1073  //cotbeta_ = gv_dot_gvy / gv_dot_gvz;
1074 }
1075 
1076 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:42
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
float xx() const
Definition: LocalError.h:24
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
unsigned int panel() const
panel id
Definition: PXFDetId.h:52
float clusterProbability(unsigned int flags=0) const
T y() const
Definition: PV2DBase.h:45
virtual bool containsBigPixelInX(const int &ixmin, const int &ixmax) const =0
T perp() const
Definition: PV3DBase.h:71
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
virtual int ncolumns() const =0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
T y() const
Definition: PV3DBase.h:62
virtual int nrows() const =0
#define abs(x)
Definition: mlp_lapack.h:159
unsigned int ladder() const
ladder id
Definition: PXBDetId.h:39
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
unsigned int blade() const
blade id
Definition: PXFDetId.h:48
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
virtual bool isItEdgePixelInX(int ixbin) const =0
float yy() const
Definition: LocalError.h:26
T sqrt(T t)
Definition: SSEVec.h:46
T z() const
Definition: PV3DBase.h:63
unsigned int module() const
det id
Definition: PXBDetId.h:43
SiPixelErrorEstimation(const edm::ParameterSet &)
unsigned int module() const
det id
Definition: PXFDetId.h:56
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual bool containsBigPixelInY(const int &iymin, const int &iymax) const =0
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
float cl
Definition: Combine.cc:71
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
Definition: DetId.h:20
void computeAnglesFromDetPosition(const SiPixelCluster &cl, const GeomDetUnit &det, float &alpha, float &beta)
tuple tracks
Definition: testEve_cfg.py:39
ClusterRef cluster() const
Definition: SiPixelRecHit.h:41
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: Handle.h:74
edm::EventID id() const
Definition: EventBase.h:56
Pixel cluster – collection of neighboring pixels above threshold.
virtual bool isItEdgePixelInY(int iybin) const =0
iterator end()
Definition: DetSetNew.h:59
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< PSimHit > associateHit(const TrackingRecHit &thit)
float y() const
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
unsigned int side() const
positive or negative id
Definition: PXFDetId.h:38
tuple cout
Definition: gather_cfg.py:121
unsigned short adc
virtual void analyze(const edm::Event &, const edm::EventSetup &)
T x() const
Definition: PV2DBase.h:44
T x() const
Definition: PV3DBase.h:61
const PositionType & position() const
std::vector< SimTrack > SimTrackContainer
float x() const
Our base class.
Definition: SiPixelRecHit.h:22
unsigned long long le
Definition: VDTMath.h:202
iterator begin()
Definition: DetSetNew.h:56