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