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
9 
10 
11 // Services
14 
15 #include "boost/multi_array.hpp"
16 
17 #include <iostream>
18 using namespace std;
19 
20 #define NEW_CPEERROR // must be constistent with base.cc, generic cc/h and genericProducer.cc
21 
22 namespace {
23  constexpr float micronsToCm = 1.0e-4;
24  const bool MYDEBUG = false;
25 }
26 
27 //-----------------------------------------------------------------------------
29 //-----------------------------------------------------------------------------
30 #ifdef NEW_CPEERROR
32  const MagneticField * mag,
33  const TrackerGeometry& geom,
34  const TrackerTopology& ttopo,
35  const SiPixelLorentzAngle * lorentzAngle,
36  const SiPixelGenErrorDBObject * genErrorDBObject,
37  const SiPixelLorentzAngle * lorentzAngleWidth=0)
38  : PixelCPEBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, 0,lorentzAngleWidth,0) {
39 #else
41  const MagneticField * mag,
42  const TrackerGeometry& geom,
43  const TrackerTopology& ttopo,
44  const SiPixelLorentzAngle * lorentzAngle,
45  const SiPixelGenErrorDBObject * genErrorDBObject,
46  const SiPixelTemplateDBObject * templateDBobject,
47  const SiPixelLorentzAngle * lorentzAngleWidth=0)
48  : PixelCPEBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, templateDBobject,lorentzAngleWidth,0) {
49 #endif
50 
51  if (theVerboseLevel > 0)
52  LogDebug("PixelCPEGeneric")
53  << " constructing a generic algorithm for ideal pixel detector.\n"
54  << " CPEGeneric:: VerboseLevel = " << theVerboseLevel;
55 
56  // Externally settable cuts
57  the_eff_charge_cut_lowX = conf.getParameter<double>("eff_charge_cut_lowX");
58  the_eff_charge_cut_lowY = conf.getParameter<double>("eff_charge_cut_lowY");
59  the_eff_charge_cut_highX = conf.getParameter<double>("eff_charge_cut_highX");
60  the_eff_charge_cut_highY = conf.getParameter<double>("eff_charge_cut_highY");
61  the_size_cutX = conf.getParameter<double>("size_cutX");
62  the_size_cutY = conf.getParameter<double>("size_cutY");
63 
64  EdgeClusterErrorX_ = conf.getParameter<double>("EdgeClusterErrorX");
65  EdgeClusterErrorY_ = conf.getParameter<double>("EdgeClusterErrorY");
66 
67  // Externally settable flags to inflate errors
68  inflate_errors = conf.getParameter<bool>("inflate_errors");
69  inflate_all_errors_no_trk_angle = conf.getParameter<bool>("inflate_all_errors_no_trk_angle");
70 
71  UseErrorsFromTemplates_ = conf.getParameter<bool>("UseErrorsFromTemplates");
72  TruncatePixelCharge_ = conf.getParameter<bool>("TruncatePixelCharge");
73  IrradiationBiasCorrection_ = conf.getParameter<bool>("IrradiationBiasCorrection");
74  DoCosmics_ = conf.getParameter<bool>("DoCosmics");
75  //LoadTemplatesFromDB_ = conf.getParameter<bool>("LoadTemplatesFromDB");
76 
77  bool isUpgrade=false;
78  if ( conf.exists("Upgrade") && conf.getParameter<bool>("Upgrade")) {
79  isUpgrade=true;
80  xerr_barrel_ln_= {0.00114,0.00104,0.00214};
81  xerr_barrel_ln_def_=0.00425;
82  yerr_barrel_ln_= {0.00299,0.00203,0.0023,0.00237,0.00233,0.00243,0.00232,0.00259,0.00176};
83  yerr_barrel_ln_def_=0.00245;
84  xerr_endcap_= {0.00151,0.000813,0.00221};
85  xerr_endcap_def_=0.00218;
86  yerr_endcap_= {0.00261,0.00107,0.00264};
87  yerr_endcap_def_=0.00357;
88 
89  if ( conf.exists("SmallPitch") && conf.getParameter<bool>("SmallPitch")) {
90  xerr_barrel_l1_= {0.00104, 0.000691, 0.00122};
91  xerr_barrel_l1_def_=0.00321;
92  yerr_barrel_l1_= {0.00199,0.00136,0.0015,0.00153,0.00152,0.00171,0.00154,0.00157,0.00154};
93  yerr_barrel_l1_def_=0.00164;
94  }
95  else{
96  xerr_barrel_l1_= {0.00114,0.00104,0.00214};
97  xerr_barrel_l1_def_=0.00425;
98  yerr_barrel_l1_= {0.00299,0.00203,0.0023,0.00237,0.00233,0.00243,0.00232,0.00259,0.00176};
99  yerr_barrel_l1_def_=0.00245;
100  }
101  }
102  isUpgrade_=isUpgrade;
103 
104  // Select the position error source
105  // For upgrde and cosmics force the use simple errors
106  if( isUpgrade_ || (DoCosmics_) ) UseErrorsFromTemplates_ = false;
107 
108  if ( !UseErrorsFromTemplates_ && ( TruncatePixelCharge_ ||
109  IrradiationBiasCorrection_ ||
110  DoCosmics_ ||
112  throw cms::Exception("PixelCPEGeneric::PixelCPEGeneric: ")
113  << "\nERROR: UseErrorsFromTemplates_ is set to False in PixelCPEGeneric_cfi.py. "
114  << " In this case it does not make sense to set any of the following to True: "
115  << " TruncatePixelCharge_, IrradiationBiasCorrection_, DoCosmics_, LoadTemplatesFromDB_ !!!"
116  << "\n\n";
117  }
118 
119  // Use errors from templates or from GenError
120  if ( UseErrorsFromTemplates_ ) {
121 
122 #ifdef NEW_CPEERROR
123 
124  if ( LoadTemplatesFromDB_ ) { // From DB
126  throw cms::Exception("InvalidCalibrationLoaded")
127  << "ERROR: GenErrors not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
128  << ( *genErrorDBObject_ ).version();
129  if(MYDEBUG) cout<<"Loaded genErrorDBObject v"<<( *genErrorDBObject_ ).version()<<endl;
130  } else { // From file
132  throw cms::Exception("InvalidCalibrationLoaded")
133  << "ERROR: GenErrors not loaded correctly from text file. Reconstruction will fail.";
134  } // if load from DB
135 
136 #else
137 
138  if ( LoadTemplatesFromDB_ ) {
139  // Initialize template store to the selected ID [Morris, 6/25/08]
140  if ( !SiPixelTemplate::pushfile( *templateDBobject_, thePixelTemp_) )
141  throw cms::Exception("InvalidCalibrationLoaded")
142  << "ERROR: Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
143  << ( *templateDBobject_ ).version();
144  if(MYDEBUG) cout<<"Loaded templateDBobject "<<( *templateDBobject_ ).version()<<endl;
145 
146  } else {
147  if ( !SiPixelTemplate::pushfile( -999, thePixelTemp_ ) )
148  throw cms::Exception("InvalidCalibrationLoaded")
149  << "ERROR: Templates not loaded correctly from text file. Reconstruction will fail.";
150  } // if load from DB
151 
152 #endif // NEW_CPEERROR
153 
154  } // if ( UseErrorsFromTemplates_ )
155 
156  if(MYDEBUG) {
157  cout << "From PixelCPEGeneric::PixelCPEGeneric(...)" << endl;
158  cout << "(int)UseErrorsFromTemplates_ = " << (int)UseErrorsFromTemplates_ << endl;
159  cout << "TruncatePixelCharge_ = " << (int)TruncatePixelCharge_ << endl;
160  cout << "IrradiationBiasCorrection_ = " << (int)IrradiationBiasCorrection_ << endl;
161  cout << "(int)DoCosmics_ = " << (int)DoCosmics_ << endl;
162  cout << "(int)LoadTemplatesFromDB_ = " << (int)LoadTemplatesFromDB_ << endl;
163  }
164 
165 
166  // Default case for rechit errors in case other, more correct, errors are not vailable
167  // This are constants. Maybe there is a more efficienct way to store them.
168  xerr_barrel_l1_= {0.00115, 0.00120, 0.00088};
169  xerr_barrel_l1_def_=0.01030;
170  yerr_barrel_l1_= {0.00375,0.00230,0.00250,0.00250,0.00230,0.00230,0.00210,0.00210,0.00240};
171  yerr_barrel_l1_def_=0.00210;
172  xerr_barrel_ln_= {0.00115, 0.00120, 0.00088};
173  xerr_barrel_ln_def_=0.01030;
174  yerr_barrel_ln_= {0.00375,0.00230,0.00250,0.00250,0.00230,0.00230,0.00210,0.00210,0.00240};
175  yerr_barrel_ln_def_=0.00210;
176  xerr_endcap_= {0.0020, 0.0020};
177  xerr_endcap_def_=0.0020;
178  yerr_endcap_= {0.00210};
179  yerr_endcap_def_=0.00075;
180 
181 }
182 
184 {
185  return new ClusterParamGeneric(cl);
186 }
187 
188 
189 
190 //-----------------------------------------------------------------------------
194 //-----------------------------------------------------------------------------
196 PixelCPEGeneric::localPosition(DetParam const & theDetParam, ClusterParam & theClusterParamBase) const
197 {
198 
199  ClusterParamGeneric & theClusterParam = static_cast<ClusterParamGeneric &>(theClusterParamBase);
200 
201  //cout<<" in PixelCPEGeneric:localPosition - "<<endl; //dk
202 
203  float chargeWidthX = (theDetParam.lorentzShiftInCmX * theDetParam.widthLAFractionX);
204  float chargeWidthY = (theDetParam.lorentzShiftInCmY * theDetParam.widthLAFractionY);
205  float shiftX = 0.5f*theDetParam.lorentzShiftInCmX;
206  float shiftY = 0.5f*theDetParam.lorentzShiftInCmY;
207 
208  //cout<<" main la width "<<chargeWidthX<<" "<<chargeWidthY<<endl;
209 
210  if ( UseErrorsFromTemplates_ ) {
211 
212  float qclus = theClusterParam.theCluster->charge();
213  float locBz = theDetParam.bz;
214  //cout << "PixelCPEGeneric::localPosition(...) : locBz = " << locBz << endl;
215 
216  theClusterParam.pixmx = -999.9; // max pixel charge for truncation of 2-D cluster
217  theClusterParam.sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
218  theClusterParam.deltay = -999.9; // CPE Generic y-bias for multi-pixel cluster
219  theClusterParam.sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
220  theClusterParam.deltax = -999.9; // CPE Generic x-bias for multi-pixel cluster
221  theClusterParam.sy1 = -999.9; // CPE Generic y-error for single single-pixel
222  theClusterParam.dy1 = -999.9; // CPE Generic y-bias for single single-pixel cluster
223  theClusterParam.sy2 = -999.9; // CPE Generic y-error for single double-pixel cluster
224  theClusterParam.dy2 = -999.9; // CPE Generic y-bias for single double-pixel cluster
225  theClusterParam.sx1 = -999.9; // CPE Generic x-error for single single-pixel cluster
226  theClusterParam.dx1 = -999.9; // CPE Generic x-bias for single single-pixel cluster
227  theClusterParam.sx2 = -999.9; // CPE Generic x-error for single double-pixel cluster
228  theClusterParam.dx2 = -999.9; // CPE Generic x-bias for single double-pixel cluster
229 
230 
231 #ifdef NEW_CPEERROR
232 
234  int gtemplID_ = theDetParam.detTemplateId;
235 
236  //int gtemplID0 = genErrorDBObject_->getGenErrorID(theDetParam.theDet->geographicalId().rawId());
237  //if(gtemplID0!=gtemplID_) cout<<" different id "<< gtemplID_<<" "<<gtemplID0<<endl;
238 
239  theClusterParam.qBin_ = gtempl.qbin( gtemplID_, theClusterParam.cotalpha, theClusterParam.cotbeta, locBz, qclus,
240  theClusterParam.pixmx, theClusterParam.sigmay, theClusterParam.deltay,
241  theClusterParam.sigmax, theClusterParam.deltax, theClusterParam.sy1,
242  theClusterParam.dy1, theClusterParam.sy2, theClusterParam.dy2, theClusterParam.sx1,
243  theClusterParam.dx1, theClusterParam.sx2, theClusterParam.dx2 );
244 
245  // now use the charge widths stored in the new generic template headers (change to the
246  // incorrect sign convention of the base class)
247  bool useLAWidthFromGenError = false;
248  if(useLAWidthFromGenError) {
249  chargeWidthX = (-micronsToCm*gtempl.lorxwidth());
250  chargeWidthY = (-micronsToCm*gtempl.lorywidth());
251  if(MYDEBUG) cout<< " redefine la width (gen-error) "<< chargeWidthX<<" "<< chargeWidthY <<endl;
252  }
253  if(MYDEBUG) cout<<" GenError: "<<gtemplID_<<endl;
254 
255 #else // select templates
256 
257  SiPixelTemplate templ(thePixelTemp_);
258  int templID_ = theDetParam.detTemplateId;
259 
260  theClusterParam.qBin_ = templ.qbin( templID_, theClusterParam.cotalpha, theClusterParam.cotbeta, locBz, qclus, // inputs
261  theClusterParam.pixmx, // returned by reference
262  theClusterParam.sigmay, theClusterParam.deltay, theClusterParam.sigmax, theClusterParam.deltax, // returned by reference
263  theClusterParam.sy1, theClusterParam.dy1, theClusterParam.sy2, theClusterParam.dy2, theClusterParam.sx1, theClusterParam.dx1, theClusterParam.sx2, theClusterParam.dx2 ); // returned by reference
264 
265  //if(MYDEBUG) {
266  //cout<<" errors "<<templID_<<" "<<theClusterParam.qBin_<<" "<<theClusterParam.sigmax<<" "<<theClusterParam.sigmay<<endl;
267  //int templID0 = templateDBobject_->getTemplateID(theDetParam.theDet->geographicalId().rawId());
268  //if(templID0!=templID_) cout<<" different id"<< templID_<<" "<<templID0<<endl;
269  //}
270 
271 #endif // NEW_CPEERROR
272 
273  // These numbers come in microns from the qbin(...) call. Transform them to cm.
274  theClusterParam.deltax = theClusterParam.deltax * micronsToCm;
275  theClusterParam.dx1 = theClusterParam.dx1 * micronsToCm;
276  theClusterParam.dx2 = theClusterParam.dx2 * micronsToCm;
277 
278  theClusterParam.deltay = theClusterParam.deltay * micronsToCm;
279  theClusterParam.dy1 = theClusterParam.dy1 * micronsToCm;
280  theClusterParam.dy2 = theClusterParam.dy2 * micronsToCm;
281 
282  theClusterParam.sigmax = theClusterParam.sigmax * micronsToCm;
283  theClusterParam.sx1 = theClusterParam.sx1 * micronsToCm;
284  theClusterParam.sx2 = theClusterParam.sx2 * micronsToCm;
285 
286  theClusterParam.sigmay = theClusterParam.sigmay * micronsToCm;
287  theClusterParam.sy1 = theClusterParam.sy1 * micronsToCm;
288  theClusterParam.sy2 = theClusterParam.sy2 * micronsToCm;
289 
290  } // if ( UseErrorsFromTemplates_ )
291 
292  float Q_f_X = 0.0;
293  float Q_l_X = 0.0;
294  float Q_f_Y = 0.0;
295  float Q_l_Y = 0.0;
296  collect_edge_charges( theClusterParam,
297  Q_f_X, Q_l_X,
298  Q_f_Y, Q_l_Y );
299 
300  //--- Find the inner widths along X and Y in one shot. We
301  //--- compute the upper right corner of the inner pixels
302  //--- (== lower left corner of upper right pixel) and
303  //--- the lower left corner of the inner pixels
304  //--- (== upper right corner of lower left pixel), and then
305  //--- subtract these two points in the formula.
306 
307  //--- Upper Right corner of Lower Left pixel -- in measurement frame
308  MeasurementPoint meas_URcorn_LLpix( theClusterParam.theCluster->minPixelRow()+1.0,
309  theClusterParam.theCluster->minPixelCol()+1.0 );
310 
311  //--- Lower Left corner of Upper Right pixel -- in measurement frame
312  MeasurementPoint meas_LLcorn_URpix( theClusterParam.theCluster->maxPixelRow(),
313  theClusterParam.theCluster->maxPixelCol() );
314 
315  //--- These two now converted into the local
316  LocalPoint local_URcorn_LLpix;
317  LocalPoint local_LLcorn_URpix;
318 
319  // PixelCPEGeneric can be used with or without track angles
320  // If PixelCPEGeneric is called with track angles, use them to correct for bows/kinks:
321  if ( theClusterParam.with_track_angle ) {
322  local_URcorn_LLpix = theDetParam.theTopol->localPosition(meas_URcorn_LLpix, theClusterParam.loc_trk_pred);
323  local_LLcorn_URpix = theDetParam.theTopol->localPosition(meas_LLcorn_URpix, theClusterParam.loc_trk_pred);
324  } else {
325  local_URcorn_LLpix = theDetParam.theTopol->localPosition(meas_URcorn_LLpix);
326  local_LLcorn_URpix = theDetParam.theTopol->localPosition(meas_LLcorn_URpix);
327  }
328 
329  #ifdef EDM_ML_DEBUG
330  if (theVerboseLevel > 20) {
331  cout
332  << "\n\t >>> theClusterParam.theCluster->x = " << theClusterParam.theCluster->x()
333  << "\n\t >>> theClusterParam.theCluster->y = " << theClusterParam.theCluster->y()
334  << "\n\t >>> cluster: minRow = " << theClusterParam.theCluster->minPixelRow()
335  << " minCol = " << theClusterParam.theCluster->minPixelCol()
336  << "\n\t >>> cluster: maxRow = " << theClusterParam.theCluster->maxPixelRow()
337  << " maxCol = " << theClusterParam.theCluster->maxPixelCol()
338  << "\n\t >>> meas: inner lower left = " << meas_URcorn_LLpix.x()
339  << "," << meas_URcorn_LLpix.y()
340  << "\n\t >>> meas: inner upper right = " << meas_LLcorn_URpix.x()
341  << "," << meas_LLcorn_URpix.y()
342  << endl;
343  }
344 #endif
345 
346  //--- &&& Note that the cuts below should not be hardcoded (like in Orca and
347  //--- &&& CPEFromDetPosition/PixelCPEInitial), but rather be
348  //--- &&& externally settable (but tracked) parameters.
349 
350  //--- Position, including the half lorentz shift
351 
352  #ifdef EDM_ML_DEBUG
353  if (theVerboseLevel > 20)
354  cout << "\t >>> Generic:: processing X" << endl;
355 #endif
356 
357  float xPos =
358  generic_position_formula( theClusterParam.theCluster->sizeX(),
359  Q_f_X, Q_l_X,
360  local_URcorn_LLpix.x(), local_LLcorn_URpix.x(),
361  chargeWidthX, // lorentz shift in cm
362  theDetParam.theThickness,
363  theClusterParam.cotalpha,
364  theDetParam.thePitchX,
365  theDetParam.theRecTopol->isItBigPixelInX( theClusterParam.theCluster->minPixelRow() ),
366  theDetParam.theRecTopol->isItBigPixelInX( theClusterParam.theCluster->maxPixelRow() ),
369  the_size_cutX); // cut for eff charge width &&&
370 
371 
372  // apply the lorentz offset correction
373  xPos = xPos + shiftX;
374 
375 #ifdef EDM_ML_DEBUG
376  if (theVerboseLevel > 20)
377  cout << "\t >>> Generic:: processing Y" << endl;
378 #endif
379 
380  float yPos =
381  generic_position_formula( theClusterParam.theCluster->sizeY(),
382  Q_f_Y, Q_l_Y,
383  local_URcorn_LLpix.y(), local_LLcorn_URpix.y(),
384  chargeWidthY, // lorentz shift in cm
385  theDetParam.theThickness,
386  theClusterParam.cotbeta,
387  theDetParam.thePitchY,
388  theDetParam.theRecTopol->isItBigPixelInY( theClusterParam.theCluster->minPixelCol() ),
389  theDetParam.theRecTopol->isItBigPixelInY( theClusterParam.theCluster->maxPixelCol() ),
392  the_size_cutY); // cut for eff charge width &&&
393 
394  // apply the lorentz offset correction
395  yPos = yPos + shiftY;
396 
397  // Apply irradiation corrections. NOT USED FOR NOW
399  if ( theClusterParam.theCluster->sizeX() == 1 ) { // size=1
400  // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
401  //float tmp1 = (0.5 * theDetParam.lorentzShiftInCmX);
402  //cout << "Apply correction correction_dx1 = " << theClusterParam.dx1 << " to xPos = " << xPos;
403  xPos = xPos - (0.5 * theDetParam.lorentzShiftInCmX);
404  // Find if pixel is double (big).
405  bool bigInX = theDetParam.theRecTopol->isItBigPixelInX( theClusterParam.theCluster->maxPixelRow() );
406  if ( !bigInX ) xPos -= theClusterParam.dx1;
407  else xPos -= theClusterParam.dx2;
408  //cout<<" to "<<xPos<<" "<<(tmp1+theClusterParam.dx1)<<endl;
409  } else { // size>1
410  //cout << "Apply correction correction_deltax = " << theClusterParam.deltax << " to xPos = " << xPos;
411  xPos -= theClusterParam.deltax;
412  //cout<<" to "<<xPos<<endl;
413  }
414 
415  if ( theClusterParam.theCluster->sizeY() == 1 ) {
416  // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
417  yPos = yPos - (0.5 * theDetParam.lorentzShiftInCmY);
418 
419  // Find if pixel is double (big).
420  bool bigInY = theDetParam.theRecTopol->isItBigPixelInY( theClusterParam.theCluster->maxPixelCol() );
421  if ( !bigInY ) yPos -= theClusterParam.dy1;
422  else yPos -= theClusterParam.dy2;
423 
424  } else {
425  //cout << "Apply correction correction_deltay = " << theClusterParam.deltay << " to yPos = " << yPos << endl;
426  yPos -= theClusterParam.deltay;
427  }
428 
429  } // if ( IrradiationBiasCorrection_ )
430 
431  //cout<<" in PixelCPEGeneric:localPosition - pos = "<<xPos<<" "<<yPos<<endl; //dk
432 
433  //--- Now put the two together
434  LocalPoint pos_in_local( xPos, yPos );
435  return pos_in_local;
436 }
437 
438 
439 
440 //-----------------------------------------------------------------------------
445 //-----------------------------------------------------------------------------
446 float
449  float Q_f,
450  float Q_l,
451  float upper_edge_first_pix,
452  float lower_edge_last_pix,
453  float lorentz_shift,
454  float theThickness, //detector thickness
455  float cot_angle,
456  float pitch,
457  bool first_is_big,
458  bool last_is_big,
459  float eff_charge_cut_low,
460  float eff_charge_cut_high,
461  float size_cut
462  ) const
463 {
464 
465  //cout<<" in PixelCPEGeneric:generic_position_formula - "<<endl; //dk
466 
467  float geom_center = 0.5f * ( upper_edge_first_pix + lower_edge_last_pix );
468 
469  //--- The case of only one pixel in this projection is separate. Note that
470  //--- here first_pix == last_pix, so the average of the two is still the
471  //--- center of the pixel.
472  if ( size == 1 ) {return geom_center;}
473 
474  //--- Width of the clusters minus the edge (first and last) pixels.
475  //--- In the note, they are denoted x_F and x_L (and y_F and y_L)
476  float W_inner = lower_edge_last_pix - upper_edge_first_pix; // in cm
477 
478  //--- Predicted charge width from geometry
479  float W_pred = theThickness * cot_angle // geometric correction (in cm)
480  - lorentz_shift; // (in cm) &&& check fpix!
481 
482  //cout<<" in PixelCPEGeneric:generic_position_formula - "<<W_inner<<" "<<W_pred<<endl; //dk
483 
484  //--- Total length of the two edge pixels (first+last)
485  float sum_of_edge = 2.0f;
486  if (first_is_big) sum_of_edge += 1.0f;
487  if (last_is_big) sum_of_edge += 1.0f;
488 
489 
490  //--- The `effective' charge width -- particle's path in first and last pixels only
491  float W_eff = std::abs( W_pred ) - W_inner;
492 
493 
494  //--- If the observed charge width is inconsistent with the expectations
495  //--- based on the track, do *not* use W_pred-W_innner. Instead, replace
496  //--- it with an *average* effective charge width, which is the average
497  //--- length of the edge pixels.
498  //
499  // bool usedEdgeAlgo = false;
500  if ( (size >= size_cut) || (
501  ( W_eff/pitch < eff_charge_cut_low ) |
502  ( W_eff/pitch > eff_charge_cut_high ) ) )
503  {
504  W_eff = pitch * 0.5f * sum_of_edge; // ave. length of edge pixels (first+last) (cm)
505  // usedEdgeAlgo = true;
506 #ifdef EDM_ML_DEBUG
507  nRecHitsUsedEdge_++;
508 #endif
509  }
510 
511 
512  //--- Finally, compute the position in this projection
513  float Qdiff = Q_l - Q_f;
514  float Qsum = Q_l + Q_f;
515 
516  //--- Temporary fix for clusters with both first and last pixel with charge = 0
517  if(Qsum==0) Qsum=1.0f;
518  //float hit_pos = geom_center + 0.5f*(Qdiff/Qsum) * W_eff + half_lorentz_shift;
519  float hit_pos = geom_center + 0.5f*(Qdiff/Qsum) * W_eff;
520 
521  //cout<<" in PixelCPEGeneric:generic_position_formula - "<<hit_pos<<" "<<lorentz_shift*0.5<<endl; //dk
522 
523  #ifdef EDM_ML_DEBUG
524  //--- Debugging output
525 #warning "Debug printouts in PixelCPEGeneric.cc has been commented because they cannot be compiled"
526  /* This part is commented because some variables used here are not defined !!
527  if (theVerboseLevel > 20) {
528  if ( theDetParam.thePart == GeomDetEnumerators::PixelBarrel || theDetParam.thePart == GeomDetEnumerators::P1PXB ) {
529  cout << "\t >>> We are in the Barrel." ;
530  } else if ( theDetParam.thePart == GeomDetEnumerators::PixelEndcap ||
531  theDetParam.thePart == GeomDetEnumerators::P1PXEC ||
532  theDetParam.thePart == GeomDetEnumerators::P2PXEC ) {
533  cout << "\t >>> We are in the Forward." ;
534  } else {
535  cout << "\t >>> We are in an unexpected subdet " << theDetParam.thePart;
536  }
537  cout
538  << "\n\t >>> cot(angle) = " << cot_angle << " pitch = " << pitch << " size = " << size
539  << "\n\t >>> upper_edge_first_pix = " << upper_edge_first_pix
540  << "\n\t >>> lower_edge_last_pix = " << lower_edge_last_pix
541  << "\n\t >>> geom_center = " << geom_center
542  << "\n\t >>> half_lorentz_shift = " << half_lorentz_shift
543  << "\n\t >>> W_inner = " << W_inner
544  << "\n\t >>> W_pred = " << W_pred
545  << "\n\t >>> W_eff(orig) = " << fabs( W_pred ) - W_inner
546  << "\n\t >>> W_eff(used) = " << W_eff
547  << "\n\t >>> sum_of_edge = " << sum_of_edge
548  << "\n\t >>> Qdiff = " << Qdiff << " Qsum = " << Qsum
549  << "\n\t >>> hit_pos = " << hit_pos
550  << "\n\t >>> RecHits: total = " << nRecHitsTotal_
551  << " used edge = " << nRecHitsUsedEdge_
552  << endl;
553  if (usedEdgeAlgo)
554  cout << "\n\t >>> Used Edge algorithm." ;
555  else
556  cout << "\n\t >>> Used angle information." ;
557  cout << endl;
558  }
559  */
560 #endif
561 
562  return hit_pos;
563 }
564 
565 
566 //-----------------------------------------------------------------------------
570 //-----------------------------------------------------------------------------
571 void
573 collect_edge_charges(ClusterParam & theClusterParamBase,
574  float & Q_f_X,
575  float & Q_l_X,
576  float & Q_f_Y,
577  float & Q_l_Y
578  ) const
579 {
580  ClusterParamGeneric & theClusterParam = static_cast<ClusterParamGeneric &>(theClusterParamBase);
581 
582  // Initialize return variables.
583  Q_f_X = Q_l_X = 0.0;
584  Q_f_Y = Q_l_Y = 0.0;
585 
586 
587  // Obtain boundaries in index units
588  int xmin = theClusterParam.theCluster->minPixelRow();
589  int xmax = theClusterParam.theCluster->maxPixelRow();
590  int ymin = theClusterParam.theCluster->minPixelCol();
591  int ymax = theClusterParam.theCluster->maxPixelCol();
592 
593 
594  // Iterate over the pixels.
595  int isize = theClusterParam.theCluster->size();
596  for (int i = 0; i != isize; ++i)
597  {
598  auto const & pixel = theClusterParam.theCluster->pixel(i);
599  // ggiurgiu@fnal.gov: add pixel charge truncation
600  float pix_adc = pixel.adc;
602  pix_adc = std::min(pix_adc, theClusterParam.pixmx );
603 
604  //
605  // X projection
606  if ( pixel.x == xmin ) Q_f_X += pix_adc;
607  if ( pixel.x == xmax ) Q_l_X += pix_adc;
608  //
609  // Y projection
610  if ( pixel.y == ymin ) Q_f_Y += pix_adc;
611  if ( pixel.y == ymax ) Q_l_Y += pix_adc;
612  }
613 
614  return;
615 }
616 
617 
618 //============== INFLATED ERROR AND ERRORS FROM DB BELOW ================
619 
620 //-------------------------------------------------------------------------
621 // Hit error in the local frame
622 //-------------------------------------------------------------------------
623 LocalError
624 PixelCPEGeneric::localError(DetParam const & theDetParam, ClusterParam & theClusterParamBase) const
625 {
626 
627  ClusterParamGeneric & theClusterParam = static_cast<ClusterParamGeneric &>(theClusterParamBase);
628 
629  const bool localPrint = false;
630  // Default errors are the maximum error used for edge clusters.
631  // These are determined by looking at residuals for edge clusters
632  float xerr = EdgeClusterErrorX_ * micronsToCm;
633  float yerr = EdgeClusterErrorY_ * micronsToCm;
634 
635 
636  // Find if cluster is at the module edge.
637  int maxPixelCol = theClusterParam.theCluster->maxPixelCol();
638  int maxPixelRow = theClusterParam.theCluster->maxPixelRow();
639  int minPixelCol = theClusterParam.theCluster->minPixelCol();
640  int minPixelRow = theClusterParam.theCluster->minPixelRow();
641 
642  bool edgex = ( theDetParam.theRecTopol->isItEdgePixelInX( minPixelRow ) ) || ( theDetParam.theRecTopol->isItEdgePixelInX( maxPixelRow ) );
643  bool edgey = ( theDetParam.theRecTopol->isItEdgePixelInY( minPixelCol ) ) || ( theDetParam.theRecTopol->isItEdgePixelInY( maxPixelCol ) );
644 
645  unsigned int sizex = theClusterParam.theCluster->sizeX();
646  unsigned int sizey = theClusterParam.theCluster->sizeY();
647  if(MYDEBUG) {
648  if( int(sizex) != (maxPixelRow - minPixelRow+1) ) cout<<" wrong x"<<endl;
649  if( int(sizey) != (maxPixelCol - minPixelCol+1) ) cout<<" wrong y"<<endl;
650  }
651 
652  // Find if cluster contains double (big) pixels.
653  bool bigInX = theDetParam.theRecTopol->containsBigPixelInX( minPixelRow, maxPixelRow );
654  bool bigInY = theDetParam.theRecTopol->containsBigPixelInY( minPixelCol, maxPixelCol );
655 
656  if(localPrint) {
657  cout<<" endge clus "<<xerr<<" "<<yerr<<endl; //dk
658  if(bigInX || bigInY) cout<<" big "<<bigInX<<" "<<bigInY<<endl;
659  if(edgex || edgey) cout<<" edge "<<edgex<<" "<<edgey<<endl;
660  cout<<" before if "<<UseErrorsFromTemplates_<<" "<<theClusterParam.qBin_<<endl;
661  if(theClusterParam.qBin_ == 0)
662  cout<<" qbin 0! "<<edgex<<" "<<edgey<<" "<<bigInX<<" "<<bigInY<<" "
663  <<sizex<<" "<<sizey<<endl;
664  }
665 
666  //if likely(UseErrorsFromTemplates_ && (qBin_!= 0) ) {
668  //
669  // Use template errors
670  //cout << "Track angles are known. We can use either errors from templates or the error parameterization from DB." << endl;
671 
672 
673  if ( !edgex ) { // Only use this for non-edge clusters
674  if ( sizex == 1 ) {
675  if ( !bigInX ) {xerr = theClusterParam.sx1;}
676  else {xerr = theClusterParam.sx2;}
677  } else {xerr = theClusterParam.sigmax;}
678  }
679 
680  if ( !edgey ) { // Only use for non-edge clusters
681  if ( sizey == 1 ) {
682  if ( !bigInY ) {yerr = theClusterParam.sy1;}
683  else {yerr = theClusterParam.sy2;}
684  } else {yerr = theClusterParam.sigmay;}
685  }
686 
687  if(localPrint) {
688  cout<<" in if "<<edgex<<" "<<edgey<<" "<<sizex<<" "<<sizey<<endl;
689  cout<<" errors "<<xerr<<" "<<yerr<<" "<<theClusterParam.sx1<<" "<<theClusterParam.sx2<<" "<<theClusterParam.sigmax<<endl; //dk
690  }
691 
692  } else { // simple errors
693 
694  // This are the simple errors, hardcoded in the code
695  cout << "Track angles are not known and we are processing cosmics." << endl;
696  //cout << "Default angle estimation which assumes track from PV (0,0,0) does not work." << endl;
697 
698  if ( theDetParam.thePart == GeomDetEnumerators::PixelBarrel || theDetParam.thePart == GeomDetEnumerators::P1PXB ) {
699 
700  DetId id = (theDetParam.theDet->geographicalId());
701  int layer=ttopo_.layer(id);
702  if ( layer==1 ) {
703  if ( !edgex ) {
704  if ( sizex<=xerr_barrel_l1_.size() ) xerr=xerr_barrel_l1_[sizex-1];
705  else xerr=xerr_barrel_l1_def_;
706  }
707 
708  if ( !edgey ) {
709  if ( sizey<=yerr_barrel_l1_.size() ) yerr=yerr_barrel_l1_[sizey-1];
710  else yerr=yerr_barrel_l1_def_;
711  }
712  } else{ // layer 2,3
713  if ( !edgex ) {
714  if ( sizex<=xerr_barrel_ln_.size() ) xerr=xerr_barrel_ln_[sizex-1];
715  else xerr=xerr_barrel_ln_def_;
716  }
717 
718  if ( !edgey ) {
719  if ( sizey<=yerr_barrel_ln_.size() ) yerr=yerr_barrel_ln_[sizey-1];
720  else yerr=yerr_barrel_ln_def_;
721  }
722  }
723 
724  } else if ( theDetParam.thePart == GeomDetEnumerators::PixelEndcap ||
725  theDetParam.thePart == GeomDetEnumerators::P1PXEC ||
726  theDetParam.thePart == GeomDetEnumerators::P2PXEC ) { // EndCap
727 
728  if ( !edgex ) {
729  if ( sizex<=xerr_endcap_.size() ) xerr=xerr_endcap_[sizex-1];
730  else xerr=xerr_endcap_def_;
731  }
732 
733  if ( !edgey ) {
734  if ( sizey<=yerr_endcap_.size() ) yerr=yerr_endcap_[sizey-1];
735  else yerr=yerr_endcap_def_;
736  }
737  } // end endcap
738 
739  if(inflate_errors) {
740  int n_bigx = 0;
741  int n_bigy = 0;
742 
743  for (int irow = 0; irow < 7; ++irow) {
744  if ( theDetParam.theRecTopol->isItBigPixelInX( irow+minPixelRow ) ) ++n_bigx;
745  }
746 
747  for (int icol = 0; icol < 21; ++icol) {
748  if ( theDetParam.theRecTopol->isItBigPixelInY( icol+minPixelCol ) ) ++n_bigy;
749  }
750 
751  xerr = (float)(sizex + n_bigx) * theDetParam.thePitchX / std::sqrt( 12.0f );
752  yerr = (float)(sizey + n_bigy) * theDetParam.thePitchY / std::sqrt( 12.0f );
753 
754  } // if(inflate_errors)
755 
756  } // end
757 
758 #ifdef EDM_ML_DEBUG
759  if ( !(xerr > 0.0) )
760  throw cms::Exception("PixelCPEGeneric::localError")
761  << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
762 
763  if ( !(yerr > 0.0) )
764  throw cms::Exception("PixelCPEGeneric::localError")
765  << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
766 #endif
767 
768  //if(localPrint) {
769  //cout<<" errors "<<xerr<<" "<<yerr<<endl; //dk
770  //if(theClusterParam.qBin_ == 0) cout<<" qbin 0 "<<xerr<<" "<<yerr<<endl;
771  //}
772 
773  auto xerr_sq = xerr*xerr;
774  auto yerr_sq = yerr*yerr;
775 
776  return LocalError( xerr_sq, 0, yerr_sq );
777 
778 }
779 
780 
781 
782 
783 
784 
#define LogDebug(id)
T getParameter(std::string const &) const
bool containsBigPixelInY(int iymin, int iymax) const
std::vector< float > xerr_barrel_l1_
int i
Definition: DBlmapReader.cc:9
float charge() const
int minPixelCol() const
bool isItEdgePixelInY(int iybin) const
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
float the_eff_charge_cut_highY
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)
std::vector< float > xerr_barrel_ln_
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
virtual bool isItBigPixelInY(const int iybin) const
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:90
bool isItEdgePixelInX(int ixbin) const
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore > &thePixelTemp_)
T y() const
Definition: PV3DBase.h:63
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< SiPixelGenErrorStore > thePixelGenError_
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:60
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:65
bool IrradiationBiasCorrection_
#define constexpr
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:63
std::vector< float > yerr_endcap_
int maxPixelRow() const
double xPos
PixelCPEGeneric(edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelGenErrorDBObject *, const SiPixelLorentzAngle *)
The constructor.
#define likely(x)
const PixelTopology * theTopol
Definition: PixelCPEBase.h:62
int minPixelRow() const
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:255
float the_eff_charge_cut_lowX
T sqrt(T t)
Definition: SSEVec.h:48
bool inflate_all_errors_no_trk_angle
void collect_edge_charges(ClusterParam &theClusterParam, float &Q_f_X, float &Q_l_X, float &Q_f_Y, float &Q_l_Y) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:251
double yPos
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:77
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
float lorxwidth()
signed lorentz x-width (microns)
float generic_position_formula(int size, float Q_f, float Q_l, float upper_edge_first_pix, float lower_edge_last_pix, float lorentz_shift, float theThickness, float cot_angle, float pitch, bool first_is_big, bool last_is_big, float eff_charge_cut_low, float eff_charge_cut_high, float size_cut) const
bool containsBigPixelInX(int ixmin, int ixmax) const
std::vector< float > xerr_endcap_
tuple conf
Definition: dbtoconf.py:185
Definition: DetId.h:18
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:103
const TrackerTopology & ttopo_
Definition: PixelCPEBase.h:243
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:248
float the_eff_charge_cut_lowY
int maxPixelCol() const
std::vector< float > yerr_barrel_l1_
unsigned int layer(const DetId &id) const
std::vector< float > yerr_barrel_ln_
int sizeY() const
Pixel cluster – collection of neighboring pixels above threshold.
LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const
Pixel pixel(int i) const
PixelCPEBase(edm::ParameterSet const &conf, const MagneticField *mag, const TrackerGeometry &geom, const TrackerTopology &ttopo, const SiPixelLorentzAngle *lorentzAngle, const SiPixelGenErrorDBObject *genErrorDBObject, const SiPixelTemplateDBObject *templateDBobject, const SiPixelLorentzAngle *lorentzAngleWidth, int flag=0)
Definition: PixelCPEBase.cc:41
float lorywidth()
signed lorentz y-width (microns)
virtual bool isItBigPixelInX(const int ixbin) const
float y() const
tuple cout
Definition: gather_cfg.py:121
static bool pushfile(int filenum, std::vector< SiPixelGenErrorStore > &thePixelTemp_)
int sizeX() const
T x() const
Definition: PV3DBase.h:62
LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) 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)
tuple size
Write out results.
float x() const
float the_eff_charge_cut_highX
ClusterParam * createClusterParam(const SiPixelCluster &cl) const
int size() const