CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelCPEGeneric.cc
Go to the documentation of this file.
2 
5 
6 // this is needed to get errors from templates
8 
9 // Services
12 
13 #include "boost/multi_array.hpp"
14 
15 #include <iostream>
16 using namespace std;
17 
18 const double HALF_PI = 1.57079632679489656;
19 
20 //-----------------------------------------------------------------------------
22 //-----------------------------------------------------------------------------
24  const MagneticField * mag, const SiPixelLorentzAngle * lorentzAngle, const SiPixelCPEGenericErrorParm * genErrorParm, const SiPixelTemplateDBObject * templateDBobject)
25  : PixelCPEBase(conf, mag, lorentzAngle, genErrorParm, templateDBobject)
26 {
27 
28  if (theVerboseLevel > 0)
29  LogDebug("PixelCPEGeneric")
30  << " constructing a generic algorithm for ideal pixel detector.\n"
31  << " CPEGeneric:: VerboseLevel = " << theVerboseLevel;
32 
33  // Externally settable cuts
34  the_eff_charge_cut_lowX = conf.getParameter<double>("eff_charge_cut_lowX");
35  the_eff_charge_cut_lowY = conf.getParameter<double>("eff_charge_cut_lowY");
36  the_eff_charge_cut_highX = conf.getParameter<double>("eff_charge_cut_highX");
37  the_eff_charge_cut_highY = conf.getParameter<double>("eff_charge_cut_highY");
38  the_size_cutX = conf.getParameter<double>("size_cutX");
39  the_size_cutY = conf.getParameter<double>("size_cutY");
40 
41  EdgeClusterErrorX_ = conf.getParameter<double>("EdgeClusterErrorX");
42  EdgeClusterErrorY_ = conf.getParameter<double>("EdgeClusterErrorY");
43 
44  // Externally settable flags to inflate errors
45  inflate_errors = conf.getParameter<bool>("inflate_errors");
46  inflate_all_errors_no_trk_angle = conf.getParameter<bool>("inflate_all_errors_no_trk_angle");
47 
48  UseErrorsFromTemplates_ = conf.getParameter<bool>("UseErrorsFromTemplates");
49  TruncatePixelCharge_ = conf.getParameter<bool>("TruncatePixelCharge");
50  IrradiationBiasCorrection_ = conf.getParameter<bool>("IrradiationBiasCorrection");
51  DoCosmics_ = conf.getParameter<bool>("DoCosmics");
52  LoadTemplatesFromDB_ = conf.getParameter<bool>("LoadTemplatesFromDB");
53 
54  if ( !UseErrorsFromTemplates_ && ( TruncatePixelCharge_ ||
55  IrradiationBiasCorrection_ ||
56  DoCosmics_ ||
57  LoadTemplatesFromDB_ ) )
58  {
59  throw cms::Exception("PixelCPEGeneric::PixelCPEGeneric: ")
60  << "\nERROR: UseErrorsFromTemplates_ is set to False in PixelCPEGeneric_cfi.py. "
61  << " In this case it does not make sense to set any of the following to True: "
62  << " TruncatePixelCharge_, IrradiationBiasCorrection_, DoCosmics_, LoadTemplatesFromDB_ !!!"
63  << "\n\n";
64  }
65 
66  if ( UseErrorsFromTemplates_ )
67  {
68  templID_ = -999;
69  if ( LoadTemplatesFromDB_ )
70  {
71  // Initialize template store to the selected ID [Morris, 6/25/08]
73  throw cms::Exception("InvalidCalibrationLoaded")
74  << "ERROR: Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
75  << ( *templateDBobject_ ).version() << ". Template ID is " << templID_;
76  }
77  else
78  {
79  if ( !templ_.pushfile( templID_ ) )
80  throw cms::Exception("InvalidCalibrationLoaded")
81  << "ERROR: Templates not loaded correctly from text file. Reconstruction will fail." << " Template ID is " << templID_;
82  }
83 
84  } // if ( UseErrorsFromTemplates_ )
85 
86  //cout << endl;
87  //cout << "From PixelCPEGeneric::PixelCPEGeneric(...)" << endl;
88  //cout << "(int)UseErrorsFromTemplates_ = " << (int)UseErrorsFromTemplates_ << endl;
89  //cout << "TruncatePixelCharge_ = " << (int)TruncatePixelCharge_ << endl;
90  //cout << "IrradiationBiasCorrection_ = " << (int)IrradiationBiasCorrection_ << endl;
91  //cout << "(int)DoCosmics_ = " << (int)DoCosmics_ << endl;
92  //cout << "(int)LoadTemplatesFromDB_ = " << (int)LoadTemplatesFromDB_ << endl;
93  //cout << endl;
94 
95 }
96 
97 
100  const GeomDetUnit & det) const
101 {
102  LocalPoint lp = localPosition(cluster,det);
103 
104  // ggiurgiu@jhu.edu 12/09/2010 : trk angles needed for bow/kink correction
105  if ( with_track_angle )
107  else
108  return theTopol->measurementPosition(lp);
109 
110 }
111 
112 
113 
114 //-----------------------------------------------------------------------------
118 //-----------------------------------------------------------------------------
121  const GeomDetUnit & det) const
122 {
123  setTheDet( det, cluster );
127  {
128  /*bool fpix; //!< barrel(false) or forward(true)
129  if ( thePart == GeomDetEnumerators::PixelBarrel )
130  fpix = false; // no, it's not forward -- it's barrel
131  else
132  fpix = true; // yes, it's forward
133  */
134  float qclus = cluster.charge();
135 
137 
138  Frame detFrame( theDet->surface().position(), theDet->surface().rotation() );
139  LocalVector Bfield = detFrame.toLocal( bfield );
140  float locBz = Bfield.z();
141  //cout << "PixelCPEGeneric::localPosition(...) : locBz = " << locBz << endl;
142 
143  pixmx = -999.9; // max pixel charge for truncation of 2-D cluster
144  sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
145  deltay = -999.9; // CPE Generic y-bias for multi-pixel cluster
146  sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
147  deltax = -999.9; // CPE Generic x-bias for multi-pixel cluster
148  sy1 = -999.9; // CPE Generic y-error for single single-pixel
149  dy1 = -999.9; // CPE Generic y-bias for single single-pixel cluster
150  sy2 = -999.9; // CPE Generic y-error for single double-pixel cluster
151  dy2 = -999.9; // CPE Generic y-bias for single double-pixel cluster
152  sx1 = -999.9; // CPE Generic x-error for single single-pixel cluster
153  dx1 = -999.9; // CPE Generic x-bias for single single-pixel cluster
154  sx2 = -999.9; // CPE Generic x-error for single double-pixel cluster
155  dx2 = -999.9; // CPE Generic x-bias for single double-pixel cluster
156 
157  qBin_ = templ_.qbin( templID_, cotalpha_, cotbeta_, locBz, qclus, // inputs
158  pixmx, // returned by reference
159  sigmay, deltay, sigmax, deltax, // returned by reference
160  sy1, dy1, sy2, dy2, sx1, dx1, sx2, dx2 ); // returned by reference
161 
162  // These numbers come in microns from the qbin(...) call. Transform them to cm.
163  const float micronsToCm = 1.0e-4;
164 
166  dx1 = dx1 * micronsToCm;
167  dx2 = dx2 * micronsToCm;
168 
170  dy1 = dy1 * micronsToCm;
171  dy2 = dy2 * micronsToCm;
172 
174  sx1 = sx1 * micronsToCm;
175  sx2 = sx2 * micronsToCm;
176 
178  sy1 = sy1 * micronsToCm;
179  sy2 = sy2 * micronsToCm;
180 
181  } // if ( UseErrorsFromTemplates_ )
182 
183 
184  float Q_f_X = 0.0;
185  float Q_l_X = 0.0;
186  float Q_m_X = 0.0;
187  float Q_f_Y = 0.0;
188  float Q_l_Y = 0.0;
189  float Q_m_Y = 0.0;
190  collect_edge_charges( cluster,
191  Q_f_X, Q_l_X, Q_m_X,
192  Q_f_Y, Q_l_Y, Q_m_Y );
193 
194  //--- Find the inner widths along X and Y in one shot. We
195  //--- compute the upper right corner of the inner pixels
196  //--- (== lower left corner of upper right pixel) and
197  //--- the lower left corner of the inner pixels
198  //--- (== upper right corner of lower left pixel), and then
199  //--- subtract these two points in the formula.
200 
201  //--- Upper Right corner of Lower Left pixel -- in measurement frame
202  MeasurementPoint meas_URcorn_LLpix( cluster.minPixelRow()+1.0,
203  cluster.minPixelCol()+1.0 );
204 
205  //--- Lower Left corner of Upper Right pixel -- in measurement frame
206  MeasurementPoint meas_LLcorn_URpix( cluster.maxPixelRow(),
207  cluster.maxPixelCol() );
208 
209  //--- These two now converted into the local
210 
211  LocalPoint local_URcorn_LLpix;
212  LocalPoint local_LLcorn_URpix;
213 
214 
215  // PixelCPEGeneric can be used with or without track angles
216  // If PixelCPEGeneric is called with track angles, use them to correct for bows/kinks:
217  if ( with_track_angle )
218  {
219  local_URcorn_LLpix = theTopol->localPosition(meas_URcorn_LLpix, loc_trk_pred_);
220  local_LLcorn_URpix = theTopol->localPosition(meas_LLcorn_URpix, loc_trk_pred_);
221  }
222  else
223  {
224  local_URcorn_LLpix = theTopol->localPosition(meas_URcorn_LLpix);
225  local_LLcorn_URpix = theTopol->localPosition(meas_LLcorn_URpix);
226  }
227 
228  if (theVerboseLevel > 20) {
229  cout
230  << "\n\t >>> cluster.x = " << cluster.x()
231  << "\n\t >>> cluster.y = " << cluster.y()
232  << "\n\t >>> cluster: minRow = " << cluster.minPixelRow()
233  << " minCol = " << cluster.minPixelCol()
234  << "\n\t >>> cluster: maxRow = " << cluster.maxPixelRow()
235  << " maxCol = " << cluster.maxPixelCol()
236  << "\n\t >>> meas: inner lower left = " << meas_URcorn_LLpix.x()
237  << "," << meas_URcorn_LLpix.y()
238  << "\n\t >>> meas: inner upper right = " << meas_LLcorn_URpix.x()
239  << "," << meas_LLcorn_URpix.y()
240  << endl;
241  }
242 
243  //--- &&& Note that the cuts below should not be hardcoded (like in Orca and
244  //--- &&& CPEFromDetPosition/PixelCPEInitial), but rather be
245  //--- &&& externally settable (but tracked) parameters.
246 
247 
248  //--- Position, including the half lorentz shift
249  if (theVerboseLevel > 20)
250  cout << "\t >>> Generic:: processing X" << endl;
251  float xPos =
252  generic_position_formula( cluster.sizeX(),
253  Q_f_X, Q_l_X,
254  local_URcorn_LLpix.x(), local_LLcorn_URpix.x(),
255  0.5*lorentzShiftInCmX_, // 0.5 * lorentz shift in
256  cotalpha_,
257  thePitchX,
258  theTopol->isItBigPixelInX( cluster.minPixelRow() ),
259  theTopol->isItBigPixelInX( cluster.maxPixelRow() ),
262  the_size_cutX); // cut for eff charge width &&&
263 
264 
265  if (theVerboseLevel > 20)
266  cout << "\t >>> Generic:: processing Y" << endl;
267  float yPos =
268  generic_position_formula( cluster.sizeY(),
269  Q_f_Y, Q_l_Y,
270  local_URcorn_LLpix.y(), local_LLcorn_URpix.y(),
271  0.5*lorentzShiftInCmY_, // 0.5 * lorentz shift in cm
272  cotbeta_,
273  thePitchY, // 0.5 * lorentz shift (may be 0)
274  theTopol->isItBigPixelInY( cluster.minPixelCol() ),
275  theTopol->isItBigPixelInY( cluster.maxPixelCol() ),
278  the_size_cutY); // cut for eff charge width &&&
279 
280  // Apply irradiation corrections.
282  {
283  if ( cluster.sizeX() == 1 )
284  {
285  // sanity chack
286  if ( cluster.maxPixelRow() != cluster.minPixelRow() )
287  throw cms::Exception("PixelCPEGeneric::localPosition")
288  << "\nERROR: cluster.maxPixelRow() != cluster.minPixelRow() although x-size = 1 !!!!! " << "\n\n";
289 
290  // Find if pixel is double (big).
291  bool bigInX = theTopol->isItBigPixelInX( cluster.maxPixelRow() );
292 
293  if ( !bigInX )
294  {
295  //cout << "Apply correction dx1 = " << dx1 << " to xPos = " << xPos << endl;
296  xPos -= dx1;
297  }
298  else
299  {
300  //cout << "Apply correction dx2 = " << dx2 << " to xPos = " << xPos << endl;
301  xPos -= dx2;
302  }
303  }
304  else if ( cluster.sizeX() > 1 )
305  {
306  //cout << "Apply correction correction_deltax = " << deltax << " to xPos = " << xPos << endl;
307  xPos -= deltax;
308  }
309  else
310  throw cms::Exception("PixelCPEGeneric::localPosition")
311  << "\nERROR: Unphysical cluster x-size = " << (int)cluster.sizeX() << "\n\n";
312 
313  if ( cluster.sizeY() == 1 )
314  {
315  // sanity chack
316  if ( cluster.maxPixelCol() != cluster.minPixelCol() )
317  throw cms::Exception("PixelCPEGeneric::localPosition")
318  << "\nERROR: cluster.maxPixelCol() != cluster.minPixelCol() although y-size = 1 !!!!! " << "\n\n";
319 
320  // Find if pixel is double (big).
321  bool bigInY = theTopol->isItBigPixelInY( cluster.maxPixelCol() );
322 
323  if ( !bigInY )
324  {
325  //cout << "Apply correction dy1 = " << dy1 << " to yPos = " << yPos << endl;
326  yPos -= dy1;
327  }
328  else
329  {
330  //cout << "Apply correction dy2 = " << dy2 << " to yPos = " << yPos << endl;
331  yPos -= dy2;
332  }
333  }
334  else if ( cluster.sizeY() > 1 )
335  {
336  //cout << "Apply correction deltay = " << deltay << " to yPos = " << yPos << endl;
337  yPos -= deltay;
338  }
339  else
340  throw cms::Exception("PixelCPEGeneric::localPosition")
341  << "\nERROR: Unphysical cluster y-size = " << (int)cluster.sizeY() << "\n\n";
342 
343  } // if ( IrradiationBiasCorrection_ )
344 
345  //--- Now put the two together
346  LocalPoint pos_in_local( xPos, yPos );
347  return pos_in_local;
348 }
349 
350 
351 
352 //-----------------------------------------------------------------------------
357 //-----------------------------------------------------------------------------
358 double
361  double Q_f,
362  double Q_l,
363  double upper_edge_first_pix,
364  double lower_edge_last_pix,
365  double half_lorentz_shift,
366  double cot_angle,
367  double pitch,
368  bool first_is_big,
369  bool last_is_big,
370  double eff_charge_cut_low,
371  double eff_charge_cut_high,
372  double size_cut
373  ) const
374 {
375  double geom_center = 0.5 * ( upper_edge_first_pix + lower_edge_last_pix );
376 
377  //--- The case of only one pixel in this projection is separate. Note that
378  //--- here first_pix == last_pix, so the average of the two is still the
379  //--- center of the pixel.
380  if ( size == 1 )
381  {
382  // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
384  return geom_center;
385  else
386  return geom_center + half_lorentz_shift;
387  }
388 
389 
390  //--- Width of the clusters minus the edge (first and last) pixels.
391  //--- In the note, they are denoted x_F and x_L (and y_F and y_L)
392  double W_inner = lower_edge_last_pix - upper_edge_first_pix; // in cm
393 
394 
395  //--- Predicted charge width from geometry
396  double W_pred =
397  theThickness * cot_angle // geometric correction (in cm)
398  - 2 * half_lorentz_shift; // (in cm) &&& check fpix!
399 
400 
401  //--- Total length of the two edge pixels (first+last)
402  double sum_of_edge = 0.0;
403  if (first_is_big) sum_of_edge += 2.0;
404  else sum_of_edge += 1.0;
405 
406  if (last_is_big) sum_of_edge += 2.0;
407  else sum_of_edge += 1.0;
408 
409 
410  //--- The `effective' charge width -- particle's path in first and last pixels only
411  double W_eff = fabs( W_pred ) - W_inner;
412 
413 
414  //--- If the observed charge width is inconsistent with the expectations
415  //--- based on the track, do *not* use W_pred-W_innner. Instead, replace
416  //--- it with an *average* effective charge width, which is the average
417  //--- length of the edge pixels.
418  //
419  bool usedEdgeAlgo = false;
420  if (( W_eff/pitch < eff_charge_cut_low ) ||
421  ( W_eff/pitch > eff_charge_cut_high ) || (size >= size_cut))
422  {
423  W_eff = pitch * 0.5 * sum_of_edge; // ave. length of edge pixels (first+last) (cm)
424  usedEdgeAlgo = true;
426  }
427 
428 
429  //--- Finally, compute the position in this projection
430  double Qdiff = Q_l - Q_f;
431  double Qsum = Q_l + Q_f;
432 
433  //--- Temporary fix for clusters with both first and last pixel with charge = 0
434  if(Qsum==0) Qsum=1.0;
435  double hit_pos = geom_center + 0.5*(Qdiff/Qsum) * W_eff + half_lorentz_shift;
436 
437  //--- Debugging output
438  if (theVerboseLevel > 20) {
440  cout << "\t >>> We are in the Barrel." ;
441  } else {
442  cout << "\t >>> We are in the Forward." ;
443  }
444  cout
445  << "\n\t >>> cot(angle) = " << cot_angle << " pitch = " << pitch << " size = " << size
446  << "\n\t >>> upper_edge_first_pix = " << upper_edge_first_pix
447  << "\n\t >>> lower_edge_last_pix = " << lower_edge_last_pix
448  << "\n\t >>> geom_center = " << geom_center
449  << "\n\t >>> half_lorentz_shift = " << half_lorentz_shift
450  << "\n\t >>> W_inner = " << W_inner
451  << "\n\t >>> W_pred = " << W_pred
452  << "\n\t >>> W_eff(orig) = " << fabs( W_pred ) - W_inner
453  << "\n\t >>> W_eff(used) = " << W_eff
454  << "\n\t >>> sum_of_edge = " << sum_of_edge
455  << "\n\t >>> Qdiff = " << Qdiff << " Qsum = " << Qsum
456  << "\n\t >>> hit_pos = " << hit_pos
457  << "\n\t >>> RecHits: total = " << nRecHitsTotal_
458  << " used edge = " << nRecHitsUsedEdge_
459  << endl;
460  if (usedEdgeAlgo)
461  cout << "\n\t >>> Used Edge algorithm." ;
462  else
463  cout << "\n\t >>> Used angle information." ;
464  cout << endl;
465  }
466 
467 
468  return hit_pos;
469 }
470 
471 
472 
473 
474 
475 //-----------------------------------------------------------------------------
479 //-----------------------------------------------------------------------------
480 void
483  float & Q_f_X,
484  float & Q_l_X,
485  float & Q_m_X,
486  float & Q_f_Y,
487  float & Q_l_Y,
488  float & Q_m_Y
489  ) const
490 {
491  // Initialize return variables.
492  Q_f_X = Q_l_X = Q_m_X = 0.0;
493  Q_f_Y = Q_l_Y = Q_m_Y = 0.0;
494 
495  // Fetch the pixels vector from the cluster.
496  const vector<SiPixelCluster::Pixel>& pixelsVec = cluster.pixels();
497 
498  // Obtain boundaries in index units
499  int xmin = cluster.minPixelRow();
500  int xmax = cluster.maxPixelRow();
501  int ymin = cluster.minPixelCol();
502  int ymax = cluster.maxPixelCol();
503 
504 // // Obtain the cluster boundaries (note: in measurement units!)
505 // float xmin = cluster.minPixelRow()+0.5;
506 // float xmax = cluster.maxPixelRow()+0.5;
507 // float ymin = cluster.minPixelCol()+0.5;
508 // float ymax = cluster.maxPixelCol()+0.5;
509 
510  // Iterate over the pixels.
511  int isize = pixelsVec.size();
512 
513  for (int i = 0; i < isize; ++i)
514  {
515  float pix_adc = -999.9;
516 
517  // ggiurgiu@fnal.gov: add pixel charge truncation
519  pix_adc = min( (float)(pixelsVec[i].adc), pixmx );
520  else
521  pix_adc = pixelsVec[i].adc;
522 
523  //
524  // X projection
525  if ( pixelsVec[i].x == xmin ) // need to match with tolerance!!! &&&
526  Q_f_X += pix_adc;
527  else if ( pixelsVec[i].x == xmax )
528  Q_l_X += pix_adc;
529  else
530  Q_m_X += pix_adc;
531  //
532  // Y projection
533  if ( pixelsVec[i].y == ymin )
534  Q_f_Y += pix_adc;
535  else if ( pixelsVec[i].y == ymax )
536  Q_l_Y += pix_adc;
537  else
538  Q_m_Y += pix_adc;
539  }
540 
541  return;
542 }
543 
544 
545 //============== INFLATED ERROR AND ERRORS FROM DB BELOW ================
546 
547 //-------------------------------------------------------------------------
548 // Hit error in the local frame
549 //-------------------------------------------------------------------------
550 LocalError
552  const GeomDetUnit & det) const
553 {
554  setTheDet( det, cluster );
555 
556  // The squared errors
557  float xerr_sq = -99999.9;
558  float yerr_sq = -99999.9;
559 
560  int sizex = cluster.sizeX();
561  int sizey = cluster.sizeY();
562 
563  // Default errors are the maximum error used for edge clusters.
564  /*
565  int row_offset = cluster.minPixelRow();
566  int col_offset = cluster.minPixelCol();
567  int n_bigInX = 0;
568  for (int irow = 0; irow < sizex; ++irow)
569  if ( theTopol->isItBigPixelInX( irow+row_offset ) )
570  ++n_bigInX;
571  int n_bigInY = 0;
572  for (int icol = 0; icol < sizey; ++icol)
573  if ( theTopol->isItBigPixelInY( icol+col_offset ) )
574  ++n_bigInX;
575  float xerr = (float)(sizex + n_bigInX) * thePitchX / sqrt(12.0);
576  float yerr = (float)(sizey + n_bigInY) * thePitchY / sqrt(12.0);
577  */
578 
579  // These are determined by looking at residuals for edge clusters
580  const float micronsToCm = 1.0e-4;
581  float xerr = EdgeClusterErrorX_ * micronsToCm;
582  float yerr = EdgeClusterErrorY_ * micronsToCm;
583 
584  // Find if cluster is at the module edge.
585  int maxPixelCol = cluster.maxPixelCol();
586  int maxPixelRow = cluster.maxPixelRow();
587  int minPixelCol = cluster.minPixelCol();
588  int minPixelRow = cluster.minPixelRow();
589  bool edgex = ( theTopol->isItEdgePixelInX( minPixelRow ) ) || ( theTopol->isItEdgePixelInX( maxPixelRow ) );
590  bool edgey = ( theTopol->isItEdgePixelInY( minPixelCol ) ) || ( theTopol->isItEdgePixelInY( maxPixelCol ) );
591 
592  // Find if cluster contains double (big) pixels.
593  bool bigInX = theTopol->containsBigPixelInX( minPixelRow, maxPixelRow );
594  bool bigInY = theTopol->containsBigPixelInY( minPixelCol, maxPixelCol );
595 
596  if ( !with_track_angle && DoCosmics_ )
597  {
598  //cout << "Track angles are not known and we are processing cosmics." << endl;
599  //cout << "Default angle estimation which assumes track from PV (0,0,0) does not work." << endl;
600  //cout << "Use an error parameterization which only depends on cluster size (by Vincenzo Chiochia)." << endl;
601 
603  {
604  if ( !edgex )
605  {
606  if ( sizex == 1 ) xerr = 0.00115; // Size = 1 -> Sigma = 11.5 um
607  else if ( sizex == 2 ) xerr = 0.00120; // Size = 2 -> Sigma = 12.0 um
608  else if ( sizex == 3 ) xerr = 0.00088; // Size = 3 -> Sigma = 8.8 um
609  else xerr = 0.01030;
610  }
611 
612  if ( !edgey )
613  {
614  if ( sizey == 1 ) yerr = 0.00375; // 37.5 um
615  else if ( sizey == 2 ) yerr = 0.00230; // 23 um
616  else if ( sizey == 3 ) yerr = 0.00250; // 25 um
617  else if ( sizey == 4 ) yerr = 0.00250; // 25 um
618  else if ( sizey == 5 ) yerr = 0.00230; // 23 um
619  else if ( sizey == 6 ) yerr = 0.00230; // 23 um
620  else if ( sizey == 7 ) yerr = 0.00210; // 21 um
621  else if ( sizey == 8 ) yerr = 0.00210; // 21 um
622  else if ( sizey == 9 ) yerr = 0.00240; // 24 um
623  else yerr = 0.00210; // 21um
624  }
625  }
626  else // EndCap
627  {
628  if ( !edgex )
629  {
630  if ( sizex == 1 ) xerr = 0.0020;
631  else if ( sizex == 2 ) xerr = 0.0020;
632  else xerr = 0.0020;
633  }
634 
635  if ( !edgey )
636  {
637  if ( sizey == 1 ) yerr = 0.00210; // 21 um
638  else yerr = 0.00075; // 7.5 um
639  }
640  }
641 
642  } // if ( !with_track_angle )
643  else
644  {
645  //cout << "Track angles are known. We can use either errors from templates or the error parameterization from DB." << endl;
646 
648  {
649  if (qBin_ == 0 && inflate_errors )
650  {
651  int n_bigx = 0;
652  int n_bigy = 0;
653 
654  int row_offset = cluster.minPixelRow();
655  int col_offset = cluster.minPixelCol();
656 
657  for (int irow = 0; irow < 7; ++irow)
658  {
659  if ( theTopol->isItBigPixelInX( irow+row_offset ) )
660  ++n_bigx;
661  }
662 
663  for (int icol = 0; icol < 21; ++icol)
664  {
665  if ( theTopol->isItBigPixelInY( icol+col_offset ) )
666  ++n_bigy;
667  }
668 
669  xerr = (float)(sizex + n_bigx) * thePitchX / sqrt( 12.0 );
670  yerr = (float)(sizey + n_bigy) * thePitchY / sqrt( 12.0 );
671 
672  } // if ( qbin == 0 && inflate_errors )
673  else
674  {
675  // Default errors
676 
677  if ( !edgex )
678  {
679  if ( sizex == 1 )
680  {
681  if ( !bigInX )
682  xerr = sx1;
683  else
684  xerr = sx2;
685  }
686  else if ( sizex > 1 )
687  xerr = sigmax;
688  else
689  throw cms::Exception("PixelCPEGeneric::localError")
690  << "\nERROR: Unphysical cluster x-size = " << sizex << "\n\n";
691  }
692 
693  if ( !edgey )
694  {
695  if ( sizey == 1 )
696  {
697  if ( !bigInY )
698  yerr = sy1;
699  else
700  yerr = sy2;
701  }
702  else if ( sizey > 1 )
703  yerr = sigmay;
704  else
705  throw cms::Exception("PixelCPEGeneric::localError")
706  << "\nERROR: Unphysical cluster y-size = " << sizex << "\n\n";
707  }
708  } // if ( qbin == 0 && inflate_errors ) else
709 
710  } //if ( UseErrorsFromTemplates_ )
711  else
712  {
713  //cout << endl << "Use errors from DB:" << endl;
714 
715  if ( edgex && edgey )
716  {
717  //--- Both axes on the edge, no point in calling PixelErrorParameterization,
718  //--- just return the max errors on both.
719  }
720  else
721  {
722  pair<float,float> errPair =
723  genErrorsFromDB_->getError( genErrorParm_, thePart, cluster.sizeX(), cluster.sizeY(),
724  alpha_, beta_, bigInX, bigInY );
725  if ( !edgex )
726  xerr = errPair.first;
727  if ( !edgey )
728  yerr = errPair.second;
729  }
730 
731  if (theVerboseLevel > 9)
732  {
733  LogDebug("PixelCPEGeneric") <<
734  " Sizex = " << cluster.sizeX() << " Sizey = " << cluster.sizeY() <<
735  " Edgex = " << edgex << " Edgey = " << edgey <<
736  " ErrX = " << xerr << " ErrY = " << yerr;
737  }
738 
739  } //if ( UseErrorsFromTemplates_ ) else
740 
741  } // if ( !with_track_angle ) else
742 
743  if ( !(xerr > 0.0) )
744  throw cms::Exception("PixelCPEGeneric::localError")
745  << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
746 
747  if ( !(yerr > 0.0) )
748  throw cms::Exception("PixelCPEGeneric::localError")
749  << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
750 
751  xerr_sq = xerr*xerr;
752  yerr_sq = yerr*yerr;
753 
754  return LocalError( xerr_sq, 0, yerr_sq );
755 
756 }
757 
758 
759 
760 
761 
762 
int adc(sample_type sample)
get the ADC sample (12 bits)
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
Topology::LocalTrackPred loc_trk_pred_
Definition: PixelCPEBase.h:270
float charge() const
int minPixelCol() const
int qbin(int id, float cotalpha, float cotbeta, float locBz, float qclus, float &pixmx, float &sigmay, float &deltay, float &sigmax, float &deltax, float &sy1, float &dy1, float &sy2, float &dy2, float &sx1, float &dx1, float &sx2, float &dx2, float &lorywidth, float &lorxwidth)
void collect_edge_charges(const SiPixelCluster &cluster, float &Q_f_X, float &Q_l_X, float &Q_m_X, float &Q_f_Y, float &Q_l_Y, float &Q_m_Y) const
const float micronsToCm
bool with_track_angle
Definition: PixelCPEBase.h:220
virtual bool containsBigPixelInX(const int &ixmin, const int &ixmax) const =0
double generic_position_formula(int size, double Q_f, double Q_l, double upper_edge_first_pix, double lower_edge_last_pix, double half_lorentz_shift, double cot_angle, double pitch, bool first_is_big, bool last_is_big, double eff_charge_cut_low, double eff_charge_cut_high, double size_cut) const
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
SiPixelCPEGenericDBErrorParametrization * genErrorsFromDB_
short getTemplateID(const uint32_t &detid) const
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float theThickness
Definition: PixelCPEBase.h:190
T y() const
Definition: PV3DBase.h:62
#define min(a, b)
Definition: mlp_lapack.h:161
std::pair< float, float > getError(const SiPixelCPEGenericErrorParm *parmErrors, GeomDetType::SubDetector pixelPart, int sizex, int sizey, float alpha, float beta, bool bigInX=false, bool bigInY=false)
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:188
bool IrradiationBiasCorrection_
double lorentzShiftInCmY_
Definition: PixelCPEBase.h:250
int maxPixelRow() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
double xPos
virtual bool isItEdgePixelInX(int ixbin) const =0
double the_eff_charge_cut_highX
int minPixelRow() const
LocalPoint localPosition(const SiPixelCluster &cluster, const GeomDetUnit &det) const
void computeLorentzShifts() const
double the_eff_charge_cut_lowY
T sqrt(T t)
Definition: SSEVec.h:46
PixelCPEGeneric(edm::ParameterSet const &conf, const MagneticField *, const SiPixelLorentzAngle *, const SiPixelCPEGenericErrorParm *, const SiPixelTemplateDBObject *)
The constructor.
bool inflate_all_errors_no_trk_angle
T z() const
Definition: PV3DBase.h:63
double lorentzShiftInCmX_
Definition: PixelCPEBase.h:249
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:263
double yPos
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
const MagneticField * magfield_
Definition: PixelCPEBase.h:257
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
virtual bool containsBigPixelInY(const int &iymin, const int &iymax) const =0
LocalError localError(const SiPixelCluster &cl, const GeomDetUnit &det) const
MeasurementPoint measurementPosition(const SiPixelCluster &, const GeomDetUnit &det) const
SiPixelTemplate templ_
tuple conf
Definition: dbtoconf.py:185
virtual bool isItBigPixelInX(const int ixbin) const =0
bool pushfile(int filenum)
LocalTrajectoryParameters loc_traj_param_
Definition: PixelCPEBase.h:272
const PixelTopology * theTopol
Definition: PixelCPEBase.h:185
int maxPixelCol() const
int sizeY() const
void setTheDet(const GeomDetUnit &det, const SiPixelCluster &cluster) const
Definition: PixelCPEBase.cc:95
Pixel cluster – collection of neighboring pixels above threshold.
virtual bool isItEdgePixelInY(int iybin) const =0
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double the_eff_charge_cut_lowX
float y() const
const RotationType & rotation() const
tuple cout
Definition: gather_cfg.py:121
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:181
Definition: DDAxes.h:10
double the_eff_charge_cut_highY
const double HALF_PI
int sizeX() const
T x() const
Definition: PV3DBase.h:61
const PositionType & position() const
tuple size
Write out results.
float x() const
const std::vector< Pixel > pixels() const
const SiPixelCPEGenericErrorParm * genErrorParm_
Definition: PixelCPEBase.h:261
int nRecHitsUsedEdge_
Definition: PixelCPEBase.h:217
virtual bool isItBigPixelInY(const int iybin) const =0