CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PixelCPETemplateReco.cc
Go to the documentation of this file.
1 // Include our own header first
3 
4 // Geometry services
7 
8 //#define DEBUG
9 
10 // MessageLogger
12 
13 // Magnetic field
15 
16 // The template header files
18 
19 // Commented for now (3/10/17) until we figure out how to resuscitate 2D template splitter
21 
23 
24 #include <vector>
25 #include "boost/multi_array.hpp"
26 
27 #include <iostream>
28 
29 using namespace SiPixelTemplateReco;
30 //using namespace SiPixelTemplateSplit;
31 using namespace std;
32 
33 namespace {
34  constexpr float micronsToCm = 1.0e-4;
35  constexpr int cluster_matrix_size_x = 13;
36  constexpr int cluster_matrix_size_y = 21;
37 } // namespace
38 
39 //-----------------------------------------------------------------------------
40 // Constructor.
41 //
42 //-----------------------------------------------------------------------------
44  const MagneticField* mag,
45  const TrackerGeometry& geom,
46  const TrackerTopology& ttopo,
47  const SiPixelLorentzAngle* lorentzAngle,
48  const SiPixelTemplateDBObject* templateDBobject)
49  : PixelCPEBase(conf, mag, geom, ttopo, lorentzAngle, nullptr, templateDBobject, nullptr, 1) {
50  //cout << endl;
51  //cout << "Constructing PixelCPETemplateReco::PixelCPETemplateReco(...)................................................." << endl;
52  //cout << endl;
53 
54  // Configurable parameters
55  //DoCosmics_ = conf.getParameter<bool>("DoCosmics"); // Not used in templates
56  //LoadTemplatesFromDB_ = conf.getParameter<bool>("LoadTemplatesFromDB"); // Moved to Base
57 
58  //cout << " PixelCPETemplateReco : (int)LoadTemplatesFromDB_ = " << (int)LoadTemplatesFromDB_ << endl;
59  //cout << "field_magnitude = " << field_magnitude << endl;
60 
61  // configuration parameter to decide between DB or text file template access
62 
64  //cout << "PixelCPETemplateReco: Loading templates from database (DB) --------- " << endl;
65 
66  // Initialize template store to the selected ID [Morris, 6/25/08]
68  throw cms::Exception("PixelCPETemplateReco")
69  << "\nERROR: Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
70  << (*templateDBobject_).version() << "\n\n";
71  } else {
72  //cout << "PixelCPETemplateReco : Loading templates for barrel and forward from ASCII files ----------" << endl;
73  barrelTemplateID_ = conf.getParameter<int>("barrelTemplateID");
74  forwardTemplateID_ = conf.getParameter<int>("forwardTemplateID");
75  templateDir_ = conf.getParameter<int>("directoryWithTemplates");
76 
78  throw cms::Exception("PixelCPETemplateReco")
79  << "\nERROR: Template ID " << barrelTemplateID_
80  << " not loaded correctly from text file. Reconstruction will fail.\n\n";
81 
82  if (!SiPixelTemplate::pushfile(forwardTemplateID_, thePixelTemp_, templateDir_))
83  throw cms::Exception("PixelCPETemplateReco")
84  << "\nERROR: Template ID " << forwardTemplateID_
85  << " not loaded correctly from text file. Reconstruction will fail.\n\n";
86  }
87 
88  speed_ = conf.getParameter<int>("speed");
89  LogDebug("PixelCPETemplateReco::PixelCPETemplateReco:") << "Template speed = " << speed_ << "\n";
90 
91  UseClusterSplitter_ = conf.getParameter<bool>("UseClusterSplitter");
92 }
93 
94 //-----------------------------------------------------------------------------
95 // Clean up.
96 //-----------------------------------------------------------------------------
98 
99 std::unique_ptr<PixelCPEBase::ClusterParam> PixelCPETemplateReco::createClusterParam(const SiPixelCluster& cl) const {
100  return std::make_unique<ClusterParamTemplate>(cl);
101 }
102 
103 //------------------------------------------------------------------
104 // Public methods mandated by the base class.
105 //------------------------------------------------------------------
106 
107 //------------------------------------------------------------------
108 // The main call to the template code.
109 //------------------------------------------------------------------
110 LocalPoint PixelCPETemplateReco::localPosition(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
111  ClusterParamTemplate& theClusterParam = static_cast<ClusterParamTemplate&>(theClusterParamBase);
112 
114  throw cms::Exception("PixelCPETemplateReco::localPosition :") << "A non-pixel detector type in here?";
115  // barrel(false) or forward(true)
116  const bool fpix = GeomDetEnumerators::isEndcap(theDetParam.thePart);
117 
118  int ID = -9999;
119  if (LoadTemplatesFromDB_) {
120  int ID0 = templateDBobject_->getTemplateID(theDetParam.theDet->geographicalId()); // just to comapre
121  ID = theDetParam.detTemplateId;
122  if (ID0 != ID)
123  edm::LogError("PixelCPETemplateReco") << " different id" << ID << " " << ID0 << endl;
124  } else { // from asci file
125  if (!fpix)
126  ID = barrelTemplateID_; // barrel
127  else
128  ID = forwardTemplateID_; // forward
129  }
130  //cout << "PixelCPETemplateReco : ID = " << ID << endl;
131 
133 
134  // Preparing to retrieve ADC counts from the SiPixeltheClusterParam.theCluster-> In the cluster,
135  // we have the following:
136  // int minPixelRow(); // Minimum pixel index in the x direction (low edge).
137  // int maxPixelRow(); // Maximum pixel index in the x direction (top edge).
138  // int minPixelCol(); // Minimum pixel index in the y direction (left edge).
139  // int maxPixelCol(); // Maximum pixel index in the y direction (right edge).
140  // So the pixels from minPixelRow() will go into clust_array_2d[0][*],
141  // and the pixels from minPixelCol() will go into clust_array_2d[*][0].
142 
143  int row_offset = theClusterParam.theCluster->minPixelRow();
144  int col_offset = theClusterParam.theCluster->minPixelCol();
145 
146  // Store the coordinates of the center of the (0,0) pixel of the array that
147  // gets passed to PixelTempReco1D
148  // Will add these values to the output of PixelTempReco1D
149  float tmp_x = float(row_offset) + 0.5f;
150  float tmp_y = float(col_offset) + 0.5f;
151 
152  // Store these offsets (to be added later) in a LocalPoint after tranforming
153  // them from measurement units (pixel units) to local coordinates (cm)
154  //
155  //
156 
157  // In case of template reco failure, these are the lorentz drift corrections
158  // to be applied
159  float lorentzshiftX = 0.5f * theDetParam.lorentzShiftInCmX;
160  float lorentzshiftY = 0.5f * theDetParam.lorentzShiftInCmY;
161 
162  // ggiurgiu@jhu.edu 12/09/2010 : update call with trk angles needed for bow/kink corrections
163  LocalPoint lp;
164 
165  if (theClusterParam.with_track_angle)
166  lp = theDetParam.theTopol->localPosition(MeasurementPoint(tmp_x, tmp_y), theClusterParam.loc_trk_pred);
167  else {
168  edm::LogError("PixelCPETemplateReco") << "@SUB = PixelCPETemplateReco::localPosition"
169  << "Should never be here. PixelCPETemplateReco should always be called with "
170  "track angles. This is a bad error !!! ";
171 
172  lp = theDetParam.theTopol->localPosition(MeasurementPoint(tmp_x, tmp_y));
173  }
174 
175  // first compute matrix size
176  int mrow = 0, mcol = 0;
177  for (int i = 0; i != theClusterParam.theCluster->size(); ++i) {
178  auto pix = theClusterParam.theCluster->pixel(i);
179  int irow = int(pix.x);
180  int icol = int(pix.y);
181  mrow = std::max(mrow, irow);
182  mcol = std::max(mcol, icol);
183  }
184  mrow -= row_offset;
185  mrow += 1;
186  mrow = std::min(mrow, cluster_matrix_size_x);
187  mcol -= col_offset;
188  mcol += 1;
189  mcol = std::min(mcol, cluster_matrix_size_y);
190  assert(mrow > 0);
191  assert(mcol > 0);
192 
193  float clustMatrix[mrow][mcol];
194  memset(clustMatrix, 0, sizeof(float) * mrow * mcol);
195 
196  // Copy clust's pixels (calibrated in electrons) into clusMatrix;
197  for (int i = 0; i != theClusterParam.theCluster->size(); ++i) {
198  auto pix = theClusterParam.theCluster->pixel(i);
199  int irow = int(pix.x) - row_offset;
200  int icol = int(pix.y) - col_offset;
201 
202  // Gavril : what do we do here if the row/column is larger than cluster_matrix_size_x/cluster_matrix_size_y ?
203  // Ignore them for the moment...
204  if ((irow < mrow) & (icol < mcol))
205  clustMatrix[irow][icol] = float(pix.adc);
206  }
207 
208  // Make and fill the bool arrays flagging double pixels
209  bool xdouble[mrow], ydouble[mcol];
210  // x directions (shorter), rows
211  for (int irow = 0; irow < mrow; ++irow)
212  xdouble[irow] = theDetParam.theRecTopol->isItBigPixelInX(irow + row_offset);
213 
214  // y directions (longer), columns
215  for (int icol = 0; icol < mcol; ++icol)
216  ydouble[icol] = theDetParam.theRecTopol->isItBigPixelInY(icol + col_offset);
217 
218  SiPixelTemplateReco::ClusMatrix clusterPayload{&clustMatrix[0][0], xdouble, ydouble, mrow, mcol};
219 
220  // Output:
221  float nonsense = -99999.9f; // nonsense init value
222  theClusterParam.templXrec_ = theClusterParam.templYrec_ = theClusterParam.templSigmaX_ =
223  theClusterParam.templSigmaY_ = nonsense;
224  // If the template recontruction fails, we want to return 1.0 for now
225  theClusterParam.templProbY_ = theClusterParam.templProbX_ = theClusterParam.templProbQ_ = 1.0f;
226  theClusterParam.templQbin_ = 0;
227  // We have a boolean denoting whether the reco failed or not
228  theClusterParam.hasFilledProb_ = false;
229 
230  float templYrec1_ = nonsense;
231  float templXrec1_ = nonsense;
232  float templYrec2_ = nonsense;
233  float templXrec2_ = nonsense;
234 
235  // ******************************************************************
236  // Do it! Use cotalpha_ and cotbeta_ calculated in PixelCPEBase
237 
238  float locBz = theDetParam.bz;
239  float locBx = theDetParam.bx;
240 
241  theClusterParam.ierr = PixelTempReco1D(ID,
242  theClusterParam.cotalpha,
243  theClusterParam.cotbeta,
244  locBz,
245  locBx,
246  clusterPayload,
247  templ,
248  theClusterParam.templYrec_,
249  theClusterParam.templSigmaY_,
250  theClusterParam.templProbY_,
251  theClusterParam.templXrec_,
252  theClusterParam.templSigmaX_,
253  theClusterParam.templProbX_,
254  theClusterParam.templQbin_,
255  speed_,
256  theClusterParam.templProbQ_);
257 
258  // ******************************************************************
259 
260  // Check exit status
261  if UNLIKELY (theClusterParam.ierr != 0) {
262  LogDebug("PixelCPETemplateReco::localPosition")
263  << "reconstruction failed with error " << theClusterParam.ierr << "\n";
264 
265  // Template reco has failed, compute position estimates based on cluster center of gravity + Lorentz drift
266  // Future improvement would be to call generic reco instead
267 
268  // ggiurgiu@jhu.edu, 21/09/2010 : trk angles needed to correct for bows/kinks
269  if (theClusterParam.with_track_angle) {
270  theClusterParam.templXrec_ =
271  theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
272  theClusterParam.templYrec_ =
273  theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
274  } else {
275  edm::LogError("PixelCPETemplateReco") << "@SUB = PixelCPETemplateReco::localPosition"
276  << "Should never be here. PixelCPETemplateReco should always be called "
277  "with track angles. This is a bad error !!! ";
278 
279  theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x()) + lorentzshiftX;
280  theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y()) + lorentzshiftY;
281  }
282  } else if UNLIKELY (UseClusterSplitter_ && theClusterParam.templQbin_ == 0) {
283  edm::LogError("PixelCPETemplateReco") << " PixelCPETemplateReco: Qbin = 0 but using cluster splitter, we should "
284  "never be here !!!!!!!!!!!!!!!!!!!!!! \n"
285  << "(int)UseClusterSplitter_ = " << (int)UseClusterSplitter_ << endl;
286 
287  //ierr =
288  //PixelTempSplit( ID, fpix, cotalpha_, cotbeta_,
289  // clust_array_2d, ydouble, xdouble,
290  // templ,
291  // templYrec1_, templYrec2_, templSigmaY_, templProbY_,
292  // templXrec1_, templXrec2_, templSigmaX_, templProbX_,
293  // templQbin_ );
294 
295  // Commented for now (3/10/17) until we figure out how to resuscitate 2D template splitter
299 
300  theClusterParam.ierr = -123;
301  /*
302  float dchisq;
303  float templProbQ_;
304  SiPixelTemplateSplit::PixelTempSplit( ID, theClusterParam.cotalpha, theClusterParam.cotbeta,
305  clust_array_2d,
306  ydouble, xdouble,
307  templ,
308  templYrec1_, templYrec2_, theClusterParam.templSigmaY_, theClusterParam.templProbY_,
309  templXrec1_, templXrec2_, theClusterParam.templSigmaX_, theClusterParam.templProbX_,
310  theClusterParam.templQbin_,
311  templProbQ_,
312  true,
313  dchisq,
314  templ2D_ );
315 
316  */
317  if (theClusterParam.ierr != 0) {
318  // Template reco has failed, compute position estimates based on cluster center of gravity + Lorentz drift
319  // Future improvement would be to call generic reco instead
320 
321  // ggiurgiu@jhu.edu, 12/09/2010 : trk angles needed to correct for bows/kinks
322  if (theClusterParam.with_track_angle) {
323  theClusterParam.templXrec_ =
324  theDetParam.theTopol->localX(theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred) + lorentzshiftX;
325  theClusterParam.templYrec_ =
326  theDetParam.theTopol->localY(theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred) + lorentzshiftY;
327  } else {
328  edm::LogError("PixelCPETemplateReco") << "@SUB = PixelCPETemplateReco::localPosition"
329  << "Should never be here. PixelCPETemplateReco should always be called "
330  "with track angles. This is a bad error !!! ";
331  theClusterParam.templXrec_ = theDetParam.theTopol->localX(theClusterParam.theCluster->x()) + lorentzshiftX;
332  theClusterParam.templYrec_ = theDetParam.theTopol->localY(theClusterParam.theCluster->y()) + lorentzshiftY;
333  }
334  } else {
335  // go from micrometer to centimeter
336  templXrec1_ *= micronsToCm;
337  templYrec1_ *= micronsToCm;
338  templXrec2_ *= micronsToCm;
339  templYrec2_ *= micronsToCm;
340 
341  // go back to the module coordinate system
342  templXrec1_ += lp.x();
343  templYrec1_ += lp.y();
344  templXrec2_ += lp.x();
345  templYrec2_ += lp.y();
346 
347  // calculate distance from each hit to the track and choose the hit closest to the track
348  float distX1 = std::abs(templXrec1_ - theClusterParam.trk_lp_x);
349  float distX2 = std::abs(templXrec2_ - theClusterParam.trk_lp_x);
350  float distY1 = std::abs(templYrec1_ - theClusterParam.trk_lp_y);
351  float distY2 = std::abs(templYrec2_ - theClusterParam.trk_lp_y);
352  theClusterParam.templXrec_ = (distX1 < distX2 ? templXrec1_ : templXrec2_);
353  theClusterParam.templYrec_ = (distY1 < distY2 ? templYrec1_ : templYrec2_);
354  }
355  } // else if ( UseClusterSplitter_ && templQbin_ == 0 )
356 
357  else // apparenly this is the good one!
358  {
359  // go from micrometer to centimeter
360  theClusterParam.templXrec_ *= micronsToCm;
361  theClusterParam.templYrec_ *= micronsToCm;
362 
363  // go back to the module coordinate system
364  theClusterParam.templXrec_ += lp.x();
365  theClusterParam.templYrec_ += lp.y();
366 
367  // Compute the Alignment Group Corrections [template ID should already be selected from call to reco procedure]
369  // Do only if the lotentzshift has meaningfull numbers
370  if (theDetParam.lorentzShiftInCmX != 0.0 || theDetParam.lorentzShiftInCmY != 0.0) {
371  // the LA width/shift returned by templates use (+)
372  // the LA width/shift produced by PixelCPEBase for positive LA is (-)
373  // correct this by iserting (-)
374  //float temp1 = -micronsToCm*templ.lorxwidth(); // old
375  //float temp2 = -micronsToCm*templ.lorywidth(); // does not incl 1/2
376  float templateLorbiasCmX = -micronsToCm * templ.lorxbias(); // new
377  float templateLorbiasCmY = -micronsToCm * templ.lorybias(); //incl. 1/2
378  // now, correctly, we can use the difference of shifts
379  //theClusterParam.templXrec_ += 0.5*(theDetParam.lorentzShiftInCmX - templateLorbiasCmX);
380  //theClusterParam.templYrec_ += 0.5*(theDetParam.lorentzShiftInCmY - templateLorbiasCmY);
381  theClusterParam.templXrec_ += (0.5 * (theDetParam.lorentzShiftInCmX) - templateLorbiasCmX);
382  theClusterParam.templYrec_ += (0.5 * (theDetParam.lorentzShiftInCmY) - templateLorbiasCmY);
383  //cout << "Templates: la lorentz offset = "
384  // <<(0.5*(theDetParam.lorentzShiftInCmX)-templateLorbiasCmX)
385  // <<" "<<templateLorbiasCmX<<" "<<templateLorbiasCmY
386  // <<" "<<temp1<<" "<<temp2
387  // <<" "<<theDetParam.lorentzShiftInCmX
388  // <<" "<<theDetParam.lorentzShiftInCmY
389  // << endl; //dk
390  } //else {cout<<" LA is 0, disable offset corrections "<<endl;} //dk
391  } //else {cout<<" Do not do LA offset correction "<<endl;} //dk
392  }
393 
394  // Save probabilities and qBin in the quantities given to us by the base class
395  // (for which there are also inline getters). &&& templProbX_ etc. should be retired...
396  theClusterParam.probabilityX_ = theClusterParam.templProbX_;
397  theClusterParam.probabilityY_ = theClusterParam.templProbY_;
398  theClusterParam.probabilityQ_ = theClusterParam.templProbQ_;
399  theClusterParam.qBin_ = theClusterParam.templQbin_;
400 
401  if (theClusterParam.ierr == 0) // always true here
402  theClusterParam.hasFilledProb_ = true;
403 
404  return LocalPoint(theClusterParam.templXrec_, theClusterParam.templYrec_);
405 }
406 
407 //------------------------------------------------------------------
408 // localError() relies on localPosition() being called FIRST!!!
409 //------------------------------------------------------------------
410 LocalError PixelCPETemplateReco::localError(DetParam const& theDetParam, ClusterParam& theClusterParamBase) const {
411  ClusterParamTemplate& theClusterParam = static_cast<ClusterParamTemplate&>(theClusterParamBase);
412 
413  //cout << endl;
414  //cout << "Set PixelCPETemplate errors .............................................." << endl;
415 
416  //cout << "CPETemplate : " << endl;
417 
418  //--- Default is the maximum error used for edge clusters.
419  //--- (never used, in fact: let comment it out, shut up the complains of the static analyzer, and save a few CPU cycles)
420  // const float sig12 = 1./sqrt(12.0);
421  // float xerr = theDetParam.thePitchX *sig12;
422  // float yerr = theDetParam.thePitchY *sig12;
423  float xerr, yerr;
424 
425  // Check if the errors were already set at the clusters splitting level
426  if (theClusterParam.theCluster->getSplitClusterErrorX() > 0.0f &&
427  theClusterParam.theCluster->getSplitClusterErrorX() < clusterSplitMaxError_ &&
428  theClusterParam.theCluster->getSplitClusterErrorY() > 0.0f &&
429  theClusterParam.theCluster->getSplitClusterErrorY() < clusterSplitMaxError_) {
430  xerr = theClusterParam.theCluster->getSplitClusterErrorX() * micronsToCm;
431  yerr = theClusterParam.theCluster->getSplitClusterErrorY() * micronsToCm;
432 
433  //cout << "Errors set at cluster splitting level : " << endl;
434  //cout << "xerr = " << xerr << endl;
435  //cout << "yerr = " << yerr << endl;
436  } else {
437  // If errors are not split at the cluster splitting level, set the errors here
438 
439  //cout << "Errors are not split at the cluster splitting level, set the errors here : " << endl;
440 
441  int maxPixelCol = theClusterParam.theCluster->maxPixelCol();
442  int maxPixelRow = theClusterParam.theCluster->maxPixelRow();
443  int minPixelCol = theClusterParam.theCluster->minPixelCol();
444  int minPixelRow = theClusterParam.theCluster->minPixelRow();
445 
446  //--- Are we near either of the edges?
447  bool edgex = (theDetParam.theRecTopol->isItEdgePixelInX(minPixelRow) ||
448  theDetParam.theRecTopol->isItEdgePixelInX(maxPixelRow));
449  bool edgey = (theDetParam.theRecTopol->isItEdgePixelInY(minPixelCol) ||
450  theDetParam.theRecTopol->isItEdgePixelInY(maxPixelCol));
451 
452  if (theClusterParam.ierr != 0) {
453  // If reconstruction fails the hit position is calculated from cluster center of gravity
454  // corrected in x by average Lorentz drift. Assign huge errors.
455  //xerr = 10.0 * (float)theClusterParam.theCluster->sizeX() * xerr;
456  //yerr = 10.0 * (float)theClusterParam.theCluster->sizeX() * yerr;
457 
459  throw cms::Exception("PixelCPETemplateReco::localPosition :") << "A non-pixel detector type in here?";
460 
461  // Assign better errors based on the residuals for failed template cases
462  if (GeomDetEnumerators::isBarrel(theDetParam.thePart)) {
463  xerr = 55.0f * micronsToCm;
464  yerr = 36.0f * micronsToCm;
465  } else {
466  xerr = 42.0f * micronsToCm;
467  yerr = 39.0f * micronsToCm;
468  }
469 
470  //cout << "xerr = " << xerr << endl;
471  //cout << "yerr = " << yerr << endl;
472 
473  //return LocalError(xerr*xerr, 0, yerr*yerr);
474  } else if (edgex || edgey) {
475  // for edge pixels assign errors according to observed residual RMS
476  if (edgex && !edgey) {
477  xerr = xEdgeXError_ * micronsToCm;
478  yerr = xEdgeYError_ * micronsToCm;
479  } else if (!edgex && edgey) {
480  xerr = yEdgeXError_ * micronsToCm;
481  yerr = yEdgeYError_ * micronsToCm;
482  } else if (edgex && edgey) {
483  xerr = bothEdgeXError_ * micronsToCm;
484  yerr = bothEdgeYError_ * micronsToCm;
485  } else {
486  throw cms::Exception(" PixelCPETemplateReco::localError: Something wrong with pixel edge flag !!!");
487  }
488 
489  //cout << "xerr = " << xerr << endl;
490  //cout << "yerr = " << yerr << endl;
491  } else {
492  // &&& need a class const
493  //const float micronsToCm = 1.0e-4;
494 
495  xerr = theClusterParam.templSigmaX_ * micronsToCm;
496  yerr = theClusterParam.templSigmaY_ * micronsToCm;
497 
498  //cout << "xerr = " << xerr << endl;
499  //cout << "yerr = " << yerr << endl;
500 
501  // &&& should also check ierr (saved as class variable) and return
502  // &&& nonsense (another class static) if the template fit failed.
503  }
504 
505  if (theVerboseLevel > 9) {
506  LogDebug("PixelCPETemplateReco") << " Sizex = " << theClusterParam.theCluster->sizeX()
507  << " Sizey = " << theClusterParam.theCluster->sizeY() << " Edgex = " << edgex
508  << " Edgey = " << edgey << " ErrX = " << xerr << " ErrY = " << yerr;
509  }
510 
511  } // else
512 
513  if (!(xerr > 0.0f))
514  throw cms::Exception("PixelCPETemplateReco::localError")
515  << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
516 
517  if (!(yerr > 0.0f))
518  throw cms::Exception("PixelCPETemplateReco::localError")
519  << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
520 
521  //cout << "Final errors set to: " << endl;
522  //cout << "xerr = " << xerr << endl;
523  //cout << "yerr = " << yerr << endl;
524  //cout << "Out of PixelCPETemplateREco..........................................................................." << endl;
525  //cout << endl;
526 
527  return LocalError(xerr * xerr, 0, yerr * yerr);
528 }
529 
531  desc.add<int>("barrelTemplateID", 0);
532  desc.add<int>("forwardTemplateID", 0);
533  desc.add<int>("directoryWithTemplates", 0);
534  desc.add<int>("speed", -2);
535  desc.add<bool>("UseClusterSplitter", false);
536 }
std::unique_ptr< ClusterParam > createClusterParam(const SiPixelCluster &cl) const override
static constexpr float xEdgeYError_
Definition: PixelCPEBase.h:247
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
short getTemplateID(const uint32_t &detid) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
uint32_t ID
Definition: Definitions.h:24
bool isBarrel(GeomDetEnumerators::SubDetector m)
T y() const
Definition: PV3DBase.h:60
int PixelTempReco1D(int id, float cotalpha, float cotbeta, float locBz, float locBx, ClusMatrix &cluster, SiPixelTemplate &templ, float &yrec, float &sigmay, float &proby, float &xrec, float &sigmax, float &probx, int &qbin, int speed, bool deadpix, std::vector< std::pair< int, int > > &zeropix, float &probQ, int &nypix, int &nxpix)
bool doLorentzFromAlignment_
Definition: PixelCPEBase.h:241
Log< level::Error, false > LogError
assert(be >=bs)
static constexpr float bothEdgeYError_
Definition: PixelCPEBase.h:253
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:47
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:52
static constexpr float yEdgeYError_
Definition: PixelCPEBase.h:250
tuple cl
Definition: haddnano.py:49
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:50
PixelCPETemplateReco(edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelTemplateDBObject *)
virtual float localX(float mpX) const =0
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
const PixelTopology * theTopol
Definition: PixelCPEBase.h:49
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:242
bool isItBigPixelInX(const int ixbin) const override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:236
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
bool isItEdgePixelInY(int iybin) const override
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isEndcap(GeomDetEnumerators::SubDetector m)
bool isItEdgePixelInX(int ixbin) const override
static constexpr float xEdgeXError_
Definition: PixelCPEBase.h:246
float lorybias()
signed lorentz y-width (microns)
virtual float localY(float mpY) const =0
static constexpr float clusterSplitMaxError_
Definition: PixelCPEBase.h:255
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore > &pixelTemp, std::string dir="CalibTracker/SiPixelESProducers/data/")
Pixel cluster – collection of neighboring pixels above threshold.
float lorxbias()
signed lorentz x-width (microns)
static constexpr float bothEdgeXError_
Definition: PixelCPEBase.h:252
bool isItBigPixelInY(const int iybin) const override
#define UNLIKELY(x)
Definition: Likely.h:21
static constexpr float yEdgeXError_
Definition: PixelCPEBase.h:249
T x() const
Definition: PV3DBase.h:59
static void fillPSetDescription(edm::ParameterSetDescription &desc)
std::vector< SiPixelTemplateStore > thePixelTemp_
bool isTrackerPixel(GeomDetEnumerators::SubDetector m)
LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const override
#define LogDebug(id)