CMS 3D CMS Logo

PixelCPEGeneric.cc
Go to the documentation of this file.
2 
6 
7 // Pixel templates contain the rec hit error parameterizaiton
9 
10 // The generic formula
12 
13 // Services
16 
17 #include "boost/multi_array.hpp"
18 
19 #include <iostream>
20 using namespace std;
21 
22 namespace {
23  constexpr float micronsToCm = 1.0e-4;
24 } // namespace
25 
26 //-----------------------------------------------------------------------------
28 //-----------------------------------------------------------------------------
30  const MagneticField* mag,
31  const TrackerGeometry& geom,
32  const TrackerTopology& ttopo,
33  const SiPixelLorentzAngle* lorentzAngle,
34  const SiPixelGenErrorDBObject* genErrorDBObject,
35  const SiPixelLorentzAngle* lorentzAngleWidth = nullptr)
36  : PixelCPEGenericBase(conf, mag, geom, ttopo, lorentzAngle, genErrorDBObject, lorentzAngleWidth) {
37  if (theVerboseLevel > 0)
38  LogDebug("PixelCPEGeneric") << " constructing a generic algorithm for ideal pixel detector.\n"
39  << " CPEGeneric:: VerboseLevel = " << theVerboseLevel;
40 
41  // Externally settable cuts
42  the_eff_charge_cut_lowX = conf.getParameter<double>("eff_charge_cut_lowX");
43  the_eff_charge_cut_lowY = conf.getParameter<double>("eff_charge_cut_lowY");
44  the_eff_charge_cut_highX = conf.getParameter<double>("eff_charge_cut_highX");
45  the_eff_charge_cut_highY = conf.getParameter<double>("eff_charge_cut_highY");
46  the_size_cutX = conf.getParameter<double>("size_cutX");
47  the_size_cutY = conf.getParameter<double>("size_cutY");
48 
49  // Externally settable flags to inflate errors
50  inflate_errors = conf.getParameter<bool>("inflate_errors");
51  inflate_all_errors_no_trk_angle = conf.getParameter<bool>("inflate_all_errors_no_trk_angle");
52 
53  NoTemplateErrorsWhenNoTrkAngles_ = conf.getParameter<bool>("NoTemplateErrorsWhenNoTrkAngles");
54  IrradiationBiasCorrection_ = conf.getParameter<bool>("IrradiationBiasCorrection");
55  DoCosmics_ = conf.getParameter<bool>("DoCosmics");
56 
57  // Upgrade means phase 2
58  isPhase2_ = conf.getParameter<bool>("Upgrade");
59 
60  // For cosmics force the use of simple errors
61  if ((DoCosmics_))
63 
65  throw cms::Exception("PixelCPEGeneric::PixelCPEGeneric: ")
66  << "\nERROR: useErrorsFromTemplates_ is set to False in PixelCPEGeneric_cfi.py. "
67  << " In this case it does not make sense to set any of the following to True: "
68  << " truncatePixelCharge_, IrradiationBiasCorrection_, DoCosmics_, LoadTemplatesFromDB_ !!!"
69  << "\n\n";
70  }
71 
72  // Use errors from templates or from GenError
74  if (LoadTemplatesFromDB_) { // From DB
76  throw cms::Exception("InvalidCalibrationLoaded")
77  << "ERROR: GenErrors not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
78  << (*genErrorDBObject_).version();
79  LogDebug("PixelCPEGeneric") << "Loaded genErrorDBObject v" << (*genErrorDBObject_).version();
80  } else { // From file
82  throw cms::Exception("InvalidCalibrationLoaded")
83  << "ERROR: GenErrors not loaded correctly from text file. Reconstruction will fail.";
84  } // if load from DB
85 
86  } else {
87 #ifdef EDM_ML_DEBUG
88  cout << " Use simple parametrised errors " << endl;
89 #endif
90  } // if ( useErrorsFromTemplates_ )
91 
92 #ifdef EDM_ML_DEBUG
93  cout << "From PixelCPEGeneric::PixelCPEGeneric(...)" << endl;
94  cout << "(int)useErrorsFromTemplates_ = " << (int)useErrorsFromTemplates_ << endl;
95  cout << "truncatePixelCharge_ = " << (int)truncatePixelCharge_ << endl;
96  cout << "IrradiationBiasCorrection_ = " << (int)IrradiationBiasCorrection_ << endl;
97  cout << "(int)DoCosmics_ = " << (int)DoCosmics_ << endl;
98  cout << "(int)LoadTemplatesFromDB_ = " << (int)LoadTemplatesFromDB_ << endl;
99 #endif
100 }
101 
102 //-----------------------------------------------------------------------------
106 //-----------------------------------------------------------------------------
107 LocalPoint PixelCPEGeneric::localPosition(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
108  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
109 
110  //cout<<" in PixelCPEGeneric:localPosition - "<<endl; //dk
111 
112  float chargeWidthX = (theDetParam.lorentzShiftInCmX * theDetParam.widthLAFractionX);
113  float chargeWidthY = (theDetParam.lorentzShiftInCmY * theDetParam.widthLAFractionY);
114  float shiftX = 0.5f * theDetParam.lorentzShiftInCmX;
115  float shiftY = 0.5f * theDetParam.lorentzShiftInCmY;
116 
117  //cout<<" main la width "<<chargeWidthX<<" "<<chargeWidthY<<endl;
118 
119  if (useErrorsFromTemplates_) {
120  float qclus = theClusterParam.theCluster->charge();
121  float locBz = theDetParam.bz;
122  float locBx = theDetParam.bx;
123  //cout << "PixelCPEGeneric::localPosition(...) : locBz = " << locBz << endl;
124 
125  theClusterParam.pixmx = -999; // max pixel charge for truncation of 2-D cluster
126  theClusterParam.sigmay = -999.9; // CPE Generic y-error for multi-pixel cluster
127  theClusterParam.deltay = -999.9; // CPE Generic y-bias for multi-pixel cluster
128  theClusterParam.sigmax = -999.9; // CPE Generic x-error for multi-pixel cluster
129  theClusterParam.deltax = -999.9; // CPE Generic x-bias for multi-pixel cluster
130  theClusterParam.sy1 = -999.9; // CPE Generic y-error for single single-pixel
131  theClusterParam.dy1 = -999.9; // CPE Generic y-bias for single single-pixel cluster
132  theClusterParam.sy2 = -999.9; // CPE Generic y-error for single double-pixel cluster
133  theClusterParam.dy2 = -999.9; // CPE Generic y-bias for single double-pixel cluster
134  theClusterParam.sx1 = -999.9; // CPE Generic x-error for single single-pixel cluster
135  theClusterParam.dx1 = -999.9; // CPE Generic x-bias for single single-pixel cluster
136  theClusterParam.sx2 = -999.9; // CPE Generic x-error for single double-pixel cluster
137  theClusterParam.dx2 = -999.9; // CPE Generic x-bias for single double-pixel cluster
138 
139  SiPixelGenError gtempl(thePixelGenError_);
140  int gtemplID_ = theDetParam.detTemplateId;
141 
142  //int gtemplID0 = genErrorDBObject_->getGenErrorID(theDetParam.theDet->geographicalId().rawId());
143  //if(gtemplID0!=gtemplID_) cout<<" different id "<< gtemplID_<<" "<<gtemplID0<<endl;
144 
145  theClusterParam.qBin_ = gtempl.qbin(gtemplID_,
146  theClusterParam.cotalpha,
147  theClusterParam.cotbeta,
148  locBz,
149  locBx,
150  qclus,
151  IrradiationBiasCorrection_,
152  theClusterParam.pixmx,
153  theClusterParam.sigmay,
154  theClusterParam.deltay,
155  theClusterParam.sigmax,
156  theClusterParam.deltax,
157  theClusterParam.sy1,
158  theClusterParam.dy1,
159  theClusterParam.sy2,
160  theClusterParam.dy2,
161  theClusterParam.sx1,
162  theClusterParam.dx1,
163  theClusterParam.sx2,
164  theClusterParam.dx2);
165 
166  // now use the charge widths stored in the new generic template headers (change to the
167  // incorrect sign convention of the base class)
168  bool useLAWidthFromGenError = false;
169  if (useLAWidthFromGenError) {
170  chargeWidthX = (-micronsToCm * gtempl.lorxwidth());
171  chargeWidthY = (-micronsToCm * gtempl.lorywidth());
172  LogDebug("PixelCPE localPosition():") << "redefine la width (gen-error)" << chargeWidthX << chargeWidthY;
173  }
174  LogDebug("PixelCPE localPosition():") << "GenError:" << gtemplID_;
175 
176  // These numbers come in microns from the qbin(...) call. Transform them to cm.
177  theClusterParam.deltax = theClusterParam.deltax * micronsToCm;
178  theClusterParam.dx1 = theClusterParam.dx1 * micronsToCm;
179  theClusterParam.dx2 = theClusterParam.dx2 * micronsToCm;
180 
181  theClusterParam.deltay = theClusterParam.deltay * micronsToCm;
182  theClusterParam.dy1 = theClusterParam.dy1 * micronsToCm;
183  theClusterParam.dy2 = theClusterParam.dy2 * micronsToCm;
184 
185  theClusterParam.sigmax = theClusterParam.sigmax * micronsToCm;
186  theClusterParam.sx1 = theClusterParam.sx1 * micronsToCm;
187  theClusterParam.sx2 = theClusterParam.sx2 * micronsToCm;
188 
189  theClusterParam.sigmay = theClusterParam.sigmay * micronsToCm;
190  theClusterParam.sy1 = theClusterParam.sy1 * micronsToCm;
191  theClusterParam.sy2 = theClusterParam.sy2 * micronsToCm;
192 
193  } // if ( useErrorsFromTemplates_ )
194  else {
195  theClusterParam.qBin_ = 0;
196  }
197 
198  int q_f_X;
199  int q_l_X;
200  int q_f_Y;
201  int q_l_Y;
202  collect_edge_charges(theClusterParam, q_f_X, q_l_X, q_f_Y, q_l_Y, useErrorsFromTemplates_ && truncatePixelCharge_);
203 
204  //--- Find the inner widths along X and Y in one shot. We
205  //--- compute the upper right corner of the inner pixels
206  //--- (== lower left corner of upper right pixel) and
207  //--- the lower left corner of the inner pixels
208  //--- (== upper right corner of lower left pixel), and then
209  //--- subtract these two points in the formula.
210 
211  //--- Upper Right corner of Lower Left pixel -- in measurement frame
212  MeasurementPoint meas_URcorn_LLpix(theClusterParam.theCluster->minPixelRow() + 1.0,
213  theClusterParam.theCluster->minPixelCol() + 1.0);
214 
215  //--- Lower Left corner of Upper Right pixel -- in measurement frame
216  MeasurementPoint meas_LLcorn_URpix(theClusterParam.theCluster->maxPixelRow(),
217  theClusterParam.theCluster->maxPixelCol());
218 
219  //--- These two now converted into the local
220  LocalPoint local_URcorn_LLpix;
221  LocalPoint local_LLcorn_URpix;
222 
223  // PixelCPEGeneric can be used with or without track angles
224  // If PixelCPEGeneric is called with track angles, use them to correct for bows/kinks:
225  if (theClusterParam.with_track_angle) {
226  local_URcorn_LLpix = theDetParam.theTopol->localPosition(meas_URcorn_LLpix, theClusterParam.loc_trk_pred);
227  local_LLcorn_URpix = theDetParam.theTopol->localPosition(meas_LLcorn_URpix, theClusterParam.loc_trk_pred);
228  } else {
229  local_URcorn_LLpix = theDetParam.theTopol->localPosition(meas_URcorn_LLpix);
230  local_LLcorn_URpix = theDetParam.theTopol->localPosition(meas_LLcorn_URpix);
231  }
232 
233 #ifdef EDM_ML_DEBUG
234  if (theVerboseLevel > 20) {
235  cout << "\n\t >>> theClusterParam.theCluster->x = " << theClusterParam.theCluster->x()
236  << "\n\t >>> theClusterParam.theCluster->y = " << theClusterParam.theCluster->y()
237  << "\n\t >>> cluster: minRow = " << theClusterParam.theCluster->minPixelRow()
238  << " minCol = " << theClusterParam.theCluster->minPixelCol()
239  << "\n\t >>> cluster: maxRow = " << theClusterParam.theCluster->maxPixelRow()
240  << " maxCol = " << theClusterParam.theCluster->maxPixelCol()
241  << "\n\t >>> meas: inner lower left = " << meas_URcorn_LLpix.x() << "," << meas_URcorn_LLpix.y()
242  << "\n\t >>> meas: inner upper right = " << meas_LLcorn_URpix.x() << "," << meas_LLcorn_URpix.y() << endl;
243  }
244 #endif
245 
246  //--- &&& Note that the cuts below should not be hardcoded (like in Orca and
247  //--- &&& CPEFromDetPosition/PixelCPEInitial), but rather be
248  //--- &&& externally settable (but tracked) parameters.
249 
250  //--- Position, including the half lorentz shift
251 
252 #ifdef EDM_ML_DEBUG
253  if (theVerboseLevel > 20)
254  cout << "\t >>> Generic:: processing X" << endl;
255 #endif
256 
258  theClusterParam.theCluster->sizeX(),
259  q_f_X,
260  q_l_X,
261  local_URcorn_LLpix.x(),
262  local_LLcorn_URpix.x(),
263  chargeWidthX, // lorentz shift in cm
264  theDetParam.theThickness,
265  theClusterParam.cotalpha,
266  theDetParam.thePitchX,
267  theDetParam.theRecTopol->isItBigPixelInX(theClusterParam.theCluster->minPixelRow()),
268  theDetParam.theRecTopol->isItBigPixelInX(theClusterParam.theCluster->maxPixelRow()),
269  the_eff_charge_cut_lowX,
270  the_eff_charge_cut_highX,
271  the_size_cutX); // cut for eff charge width &&&
272 
273  // apply the lorentz offset correction
274  xPos = xPos + shiftX;
275 
276 #ifdef EDM_ML_DEBUG
277  if (theVerboseLevel > 20)
278  cout << "\t >>> Generic:: processing Y" << endl;
279 #endif
280 
282  theClusterParam.theCluster->sizeY(),
283  q_f_Y,
284  q_l_Y,
285  local_URcorn_LLpix.y(),
286  local_LLcorn_URpix.y(),
287  chargeWidthY, // lorentz shift in cm
288  theDetParam.theThickness,
289  theClusterParam.cotbeta,
290  theDetParam.thePitchY,
291  theDetParam.theRecTopol->isItBigPixelInY(theClusterParam.theCluster->minPixelCol()),
292  theDetParam.theRecTopol->isItBigPixelInY(theClusterParam.theCluster->maxPixelCol()),
293  the_eff_charge_cut_lowY,
294  the_eff_charge_cut_highY,
295  the_size_cutY); // cut for eff charge width &&&
296 
297  // apply the lorentz offset correction
298  yPos = yPos + shiftY;
299 
300  // Apply irradiation corrections
301  if (IrradiationBiasCorrection_) {
302  if (theClusterParam.theCluster->sizeX() == 1) { // size=1
303  // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
304  //float tmp1 = (0.5 * theDetParam.lorentzShiftInCmX);
305  //cout << "Apply correction correction_dx1 = " << theClusterParam.dx1 << " to xPos = " << xPos;
306  xPos = xPos - (0.5f * theDetParam.lorentzShiftInCmX);
307  // Find if pixel is double (big).
308  bool bigInX = theDetParam.theRecTopol->isItBigPixelInX(theClusterParam.theCluster->maxPixelRow());
309  if (!bigInX)
310  xPos -= theClusterParam.dx1;
311  else
312  xPos -= theClusterParam.dx2;
313  //cout<<" to "<<xPos<<" "<<(tmp1+theClusterParam.dx1)<<endl;
314  } else { // size>1
315  //cout << "Apply correction correction_deltax = " << theClusterParam.deltax << " to xPos = " << xPos;
316  xPos -= theClusterParam.deltax;
317  //cout<<" to "<<xPos<<endl;
318  }
319 
320  if (theClusterParam.theCluster->sizeY() == 1) {
321  // ggiurgiu@jhu.edu, 02/03/09 : for size = 1, the Lorentz shift is already accounted by the irradiation correction
322  yPos = yPos - (0.5f * theDetParam.lorentzShiftInCmY);
323 
324  // Find if pixel is double (big).
325  bool bigInY = theDetParam.theRecTopol->isItBigPixelInY(theClusterParam.theCluster->maxPixelCol());
326  if (!bigInY)
327  yPos -= theClusterParam.dy1;
328  else
329  yPos -= theClusterParam.dy2;
330 
331  } else {
332  //cout << "Apply correction correction_deltay = " << theClusterParam.deltay << " to yPos = " << yPos << endl;
333  yPos -= theClusterParam.deltay;
334  }
335 
336  } // if ( IrradiationBiasCorrection_ )
337 
338  //cout<<" in PixelCPEGeneric:localPosition - pos = "<<xPos<<" "<<yPos<<endl; //dk
339 
340  //--- Now put the two together
341  LocalPoint pos_in_local(xPos, yPos);
342  return pos_in_local;
343 }
344 
345 //============== INFLATED ERROR AND ERRORS FROM DB BELOW ================
346 
347 //-------------------------------------------------------------------------
348 // Hit error in the local frame
349 //-------------------------------------------------------------------------
350 LocalError PixelCPEGeneric::localError(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
351  ClusterParamGeneric& theClusterParam = static_cast<ClusterParamGeneric&>(theClusterParamBase);
352 
353  // local variables
354  float xerr, yerr;
355  bool edgex, edgey, bigInX, bigInY;
356  int maxPixelCol, maxPixelRow, minPixelCol, minPixelRow;
357  uint sizex, sizey;
358 
359  initializeLocalErrorVariables(xerr,
360  yerr,
361  edgex,
362  edgey,
363  bigInX,
364  bigInY,
365  maxPixelCol,
366  maxPixelRow,
367  minPixelCol,
368  minPixelRow,
369  sizex,
370  sizey,
371  theDetParam,
372  theClusterParam);
373 
374  bool useTempErrors =
375  useErrorsFromTemplates_ && (!NoTemplateErrorsWhenNoTrkAngles_ || theClusterParam.with_track_angle);
376 
377  if (int(sizex) != (maxPixelRow - minPixelRow + 1))
378  LogDebug("PixelCPEGeneric") << " wrong x";
379  if (int(sizey) != (maxPixelCol - minPixelCol + 1))
380  LogDebug("PixelCPEGeneric") << " wrong y";
381 
382  LogDebug("PixelCPEGeneric") << " edge clus " << xerr << " " << yerr; //dk
383  if (bigInX || bigInY)
384  LogDebug("PixelCPEGeneric") << " big " << bigInX << " " << bigInY;
385  if (edgex || edgey)
386  LogDebug("PixelCPEGeneric") << " edge " << edgex << " " << edgey;
387  LogDebug("PixelCPEGeneric") << " before if " << useErrorsFromTemplates_ << " " << theClusterParam.qBin_;
388  if (theClusterParam.qBin_ == 0)
389  LogDebug("PixelCPEGeneric") << " qbin 0! " << edgex << " " << edgey << " " << bigInX << " " << bigInY << " "
390  << sizex << " " << sizey;
391 
392  // from PixelCPEGenericBase
393  setXYErrors(xerr, yerr, edgex, edgey, sizex, sizey, bigInX, bigInY, useTempErrors, theDetParam, theClusterParam);
394 
395  if (!useTempErrors) {
396  LogDebug("PixelCPEGeneric") << "Track angles are not known.\n"
397  << "Default angle estimation which assumes track from PV (0,0,0) does not work.";
398  }
399 
400  if (!useTempErrors && inflate_errors) {
401  int n_bigx = 0;
402  int n_bigy = 0;
403 
404  for (int irow = 0; irow < 7; ++irow) {
405  if (theDetParam.theRecTopol->isItBigPixelInX(irow + minPixelRow))
406  ++n_bigx;
407  }
408 
409  for (int icol = 0; icol < 21; ++icol) {
410  if (theDetParam.theRecTopol->isItBigPixelInY(icol + minPixelCol))
411  ++n_bigy;
412  }
413 
414  xerr = (float)(sizex + n_bigx) * theDetParam.thePitchX / std::sqrt(12.0f);
415  yerr = (float)(sizey + n_bigy) * theDetParam.thePitchY / std::sqrt(12.0f);
416  }
417 
418 #ifdef EDM_ML_DEBUG
419  if (!(xerr > 0.0))
420  throw cms::Exception("PixelCPEGeneric::localError") << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
421 
422  if (!(yerr > 0.0))
423  throw cms::Exception("PixelCPEGeneric::localError") << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
424 #endif
425 
426  LogDebug("PixelCPEGeneric") << " errors " << xerr << " " << yerr; //dk
427  if (theClusterParam.qBin_ == 0)
428  LogDebug("PixelCPEGeneric") << " qbin 0 " << xerr << " " << yerr;
429 
430  auto xerr_sq = xerr * xerr;
431  auto yerr_sq = yerr * yerr;
432 
433  return LocalError(xerr_sq, 0, yerr_sq);
434 }
435 
438  desc.add<double>("eff_charge_cut_highX", 1.0);
439  desc.add<double>("eff_charge_cut_highY", 1.0);
440  desc.add<double>("eff_charge_cut_lowX", 0.0);
441  desc.add<double>("eff_charge_cut_lowY", 0.0);
442  desc.add<double>("size_cutX", 3.0);
443  desc.add<double>("size_cutY", 3.0);
444  desc.add<double>("EdgeClusterErrorX", 50.0);
445  desc.add<double>("EdgeClusterErrorY", 85.0);
446  desc.add<bool>("inflate_errors", false);
447  desc.add<bool>("inflate_all_errors_no_trk_angle", false);
448  desc.add<bool>("NoTemplateErrorsWhenNoTrkAngles", false);
449  desc.add<bool>("UseErrorsFromTemplates", true);
450  desc.add<bool>("TruncatePixelCharge", true);
451  desc.add<bool>("IrradiationBiasCorrection", false);
452  desc.add<bool>("DoCosmics", false);
453  desc.add<bool>("Upgrade", false);
454  desc.add<bool>("SmallPitch", false);
455 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
static bool pushfile(int filenum, std::vector< SiPixelGenErrorStore > &pixelTemp, std::string dir="")
float the_eff_charge_cut_highY
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:75
LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
int maxPixelRow() const
int sizeY() const
std::vector< SiPixelGenErrorStore > thePixelGenError_
T x() const
Definition: PV2DBase.h:43
bool IrradiationBiasCorrection_
bool NoTemplateErrorsWhenNoTrkAngles_
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:50
int minPixelRow() const
T y() const
Definition: PV2DBase.h:44
int charge() const
PixelCPEGeneric(edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelGenErrorDBObject *, const SiPixelLorentzAngle *)
The constructor.
int qbin(int id, float cotalpha, float cotbeta, float locBz, float locBx, float qclus, bool irradiationCorrections, int &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)
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const PixelTopology * theTopol
Definition: PixelCPEBase.h:49
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:242
int minPixelCol() const
bool isItBigPixelInX(const int ixbin) const override
float the_eff_charge_cut_lowX
T sqrt(T t)
Definition: SSEVec.h:19
bool inflate_all_errors_no_trk_angle
int maxPixelCol() const
int sizeX() const
double f[11][100]
float lorxwidth()
signed lorentz x-width (microns)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:89
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
static const int theVerboseLevel
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:235
float the_eff_charge_cut_lowY
float y() const
float lorywidth()
signed lorentz y-width (microns)
float x() const
bool isItBigPixelInY(const int iybin) const override
static void fillPSetDescription(edm::ParameterSetDescription &desc)
float generic_position_formula(int size, int q_f, int 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)
Definition: SiPixelUtils.cc:16
LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
#define LogDebug(id)
float the_eff_charge_cut_highX