test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
PixelCPETemplateReco Class Reference

#include <PixelCPETemplateReco.h>

Inheritance diagram for PixelCPETemplateReco:
PixelCPEBase PixelClusterParameterEstimator

Classes

struct  ClusterParamTemplate
 

Public Member Functions

 PixelCPETemplateReco (edm::ParameterSet const &conf, const MagneticField *, const TrackerGeometry &, const TrackerTopology &, const SiPixelLorentzAngle *, const SiPixelTemplateDBObject *)
 
 ~PixelCPETemplateReco ()
 
- Public Member Functions inherited from PixelCPEBase
ReturnType getParameters (const SiPixelCluster &cl, const GeomDetUnit &det) const
 
ReturnType getParameters (const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) 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)
 
- Public Member Functions inherited from PixelClusterParameterEstimator
unsigned int clusterProbComputationFlag () const
 
virtual ReturnType getParameters (const SiPixelCluster &cl, const GeomDetUnit &det, const TrajectoryStateOnSurface &tsos) const
 
virtual VLocalValues localParametersV (const SiPixelCluster &cluster, const GeomDetUnit &gd) const
 
virtual VLocalValues localParametersV (const SiPixelCluster &cluster, const GeomDetUnit &gd, TrajectoryStateOnSurface &tsos) const
 
 PixelClusterParameterEstimator ()
 
virtual ~PixelClusterParameterEstimator ()
 

Private Member Functions

ClusterParamcreateClusterParam (const SiPixelCluster &cl) const
 
LocalError localError (DetParam const &theDetParam, ClusterParam &theClusterParam) const
 
LocalPoint localPosition (DetParam const &theDetParam, ClusterParam &theClusterParam) const
 

Private Attributes

int speed_
 
std::vector< SiPixelTemplateStorethePixelTemp_
 
bool UseClusterSplitter_
 

Additional Inherited Members

- Public Types inherited from PixelClusterParameterEstimator
typedef std::pair< LocalPoint,
LocalError
LocalValues
 
using ReturnType = std::tuple< LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType >
 
typedef std::vector< LocalValuesVLocalValues
 
- Protected Types inherited from PixelCPEBase
typedef GloballyPositioned
< double > 
Frame
 
- Protected Attributes inherited from PixelCPEBase
bool alpha2Order
 
bool DoLorentz_
 
const SiPixelGenErrorDBObjectgenErrorDBObject_
 
const TrackerGeometrygeom_
 
float lAOffset_
 
float lAWidthBPix_
 
float lAWidthFPix_
 
bool LoadTemplatesFromDB_
 
const SiPixelLorentzAnglelorentzAngle_
 
const SiPixelLorentzAnglelorentzAngleWidth_
 
const MagneticFieldmagfield_
 
const SiPixelTemplateDBObjecttemplateDBobject_
 
int theFlag_
 
int theVerboseLevel
 
const TrackerTopologyttopo_
 
bool useLAOffsetFromConfig_
 
bool useLAWidthFromConfig_
 
bool useLAWidthFromDB_
 
- Protected Attributes inherited from PixelClusterParameterEstimator
unsigned int clusterProbComputationFlag_
 

Detailed Description

Definition at line 34 of file PixelCPETemplateReco.h.

Constructor & Destructor Documentation

PixelCPETemplateReco::PixelCPETemplateReco ( edm::ParameterSet const &  conf,
const MagneticField mag,
const TrackerGeometry geom,
const TrackerTopology ttopo,
const SiPixelLorentzAngle lorentzAngle,
const SiPixelTemplateDBObject templateDBobject 
)

Definition at line 42 of file PixelCPETemplateReco.cc.

References Exception, edm::ParameterSet::getParameter(), PixelCPEBase::LoadTemplatesFromDB_, LogDebug, SiPixelTemplate::pushfile(), speed_, PixelCPEBase::templateDBobject_, thePixelTemp_, and UseClusterSplitter_.

48  : PixelCPEBase(conf, mag, geom, ttopo, lorentzAngle, 0, templateDBobject, 0,1)
49 {
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  {
65  //cout << "PixelCPETemplateReco: Loading templates from database (DB) --------- " << endl;
66 
67  // Initialize template store to the selected ID [Morris, 6/25/08]
69  throw cms::Exception("PixelCPETemplateReco")
70  << "\nERROR: Templates not filled correctly. Check the sqlite file. Using SiPixelTemplateDBObject version "
71  << (*templateDBobject_).version() << "\n\n";
72  }
73  else
74  {
75  //cout << "PixelCPETemplateReco : Loading templates 40 and 41 from ASCII files ------------------------" << endl;
76 
78  throw cms::Exception("PixelCPETemplateReco")
79  << "\nERROR: Templates 40 not loaded correctly from text file. Reconstruction will fail.\n\n";
80 
82  throw cms::Exception("PixelCPETemplateReco")
83  << "\nERROR: Templates 41 not loaded correctly from text file. Reconstruction will fail.\n\n";
84  }
85 
86  speed_ = conf.getParameter<int>( "speed");
87  LogDebug("PixelCPETemplateReco::PixelCPETemplateReco:") <<
88  "Template speed = " << speed_ << "\n";
89 
90  UseClusterSplitter_ = conf.getParameter<bool>("UseClusterSplitter");
91 
92 }
#define LogDebug(id)
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore > &thePixelTemp_)
std::vector< SiPixelTemplateStore > thePixelTemp_
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:255
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:251
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
PixelCPETemplateReco::~PixelCPETemplateReco ( )

Definition at line 97 of file PixelCPETemplateReco.cc.

98 {
99  // &&& delete template store?
100 }

Member Function Documentation

PixelCPEBase::ClusterParam * PixelCPETemplateReco::createClusterParam ( const SiPixelCluster cl) const
privatevirtual

Implements PixelCPEBase.

Definition at line 102 of file PixelCPETemplateReco.cc.

103 {
104  return new ClusterParamTemplate(cl);
105 }
LocalError PixelCPETemplateReco::localError ( DetParam const &  theDetParam,
ClusterParam theClusterParam 
) const
privatevirtual

Implements PixelCPEBase.

Definition at line 457 of file PixelCPETemplateReco.cc.

References Exception, f, PixelCPETemplateReco::ClusterParamTemplate::ierr, GeomDetEnumerators::isBarrel(), RectangularPixelTopology::isItEdgePixelInX(), RectangularPixelTopology::isItEdgePixelInY(), GeomDetEnumerators::isTrackerPixel(), LogDebug, mathSSE::sqrt(), PixelCPETemplateReco::ClusterParamTemplate::templSigmaX_, PixelCPETemplateReco::ClusterParamTemplate::templSigmaY_, PixelCPEBase::DetParam::thePart, PixelCPEBase::DetParam::thePitchX, PixelCPEBase::DetParam::thePitchY, PixelCPEBase::DetParam::theRecTopol, and PixelCPEBase::theVerboseLevel.

458 {
459 
460  ClusterParamTemplate & theClusterParam = static_cast<ClusterParamTemplate &>(theClusterParamBase);
461 
462  //cout << endl;
463  //cout << "Set PixelCPETemplate errors .............................................." << endl;
464 
465  //cout << "CPETemplate : " << endl;
466 
467  //--- Default is the maximum error used for edge clusters.
468  const float sig12 = 1./sqrt(12.0);
469  float xerr = theDetParam.thePitchX *sig12;
470  float yerr = theDetParam.thePitchY *sig12;
471 
472  // Check if the errors were already set at the clusters splitting level
473  if ( theClusterParam.theCluster->getSplitClusterErrorX() > 0.0f && theClusterParam.theCluster->getSplitClusterErrorX() < 7777.7f &&
474  theClusterParam.theCluster->getSplitClusterErrorY() > 0.0f && theClusterParam.theCluster->getSplitClusterErrorY() < 7777.7f )
475  {
476  xerr = theClusterParam.theCluster->getSplitClusterErrorX() * micronsToCm;
477  yerr = theClusterParam.theCluster->getSplitClusterErrorY() * micronsToCm;
478 
479  //cout << "Errors set at cluster splitting level : " << endl;
480  //cout << "xerr = " << xerr << endl;
481  //cout << "yerr = " << yerr << endl;
482  }
483  else
484  {
485  // If errors are not split at the cluster splitting level, set the errors here
486 
487  //cout << "Errors are not split at the cluster splitting level, set the errors here : " << endl;
488 
489  int maxPixelCol = theClusterParam.theCluster->maxPixelCol();
490  int maxPixelRow = theClusterParam.theCluster->maxPixelRow();
491  int minPixelCol = theClusterParam.theCluster->minPixelCol();
492  int minPixelRow = theClusterParam.theCluster->minPixelRow();
493 
494  //--- Are we near either of the edges?
495  bool edgex = ( theDetParam.theRecTopol->isItEdgePixelInX( minPixelRow ) || theDetParam.theRecTopol->isItEdgePixelInX( maxPixelRow ) );
496  bool edgey = ( theDetParam.theRecTopol->isItEdgePixelInY( minPixelCol ) || theDetParam.theRecTopol->isItEdgePixelInY( maxPixelCol ) );
497 
498  if ( theClusterParam.ierr !=0 )
499  {
500  // If reconstruction fails the hit position is calculated from cluster center of gravity
501  // corrected in x by average Lorentz drift. Assign huge errors.
502  //xerr = 10.0 * (float)theClusterParam.theCluster->sizeX() * xerr;
503  //yerr = 10.0 * (float)theClusterParam.theCluster->sizeX() * yerr;
504 
505  if(!GeomDetEnumerators::isTrackerPixel(theDetParam.thePart))
506  throw cms::Exception("PixelCPETemplateReco::localPosition :")
507  << "A non-pixel detector type in here?";
508 
509  // Assign better errors based on the residuals for failed template cases
510  if ( GeomDetEnumerators::isBarrel(theDetParam.thePart) )
511  {
512  xerr = 55.0f * micronsToCm;
513  yerr = 36.0f * micronsToCm;
514  }
515  else
516  {
517  xerr = 42.0f * micronsToCm;
518  yerr = 39.0f * micronsToCm;
519  }
520 
521  //cout << "xerr = " << xerr << endl;
522  //cout << "yerr = " << yerr << endl;
523 
524  //return LocalError(xerr*xerr, 0, yerr*yerr);
525  }
526  else if ( edgex || edgey )
527  {
528  // for edge pixels assign errors according to observed residual RMS
529  if ( edgex && !edgey )
530  {
531  xerr = 23.0f * micronsToCm;
532  yerr = 39.0f * micronsToCm;
533  }
534  else if ( !edgex && edgey )
535  {
536  xerr = 24.0f * micronsToCm;
537  yerr = 96.0f * micronsToCm;
538  }
539  else if ( edgex && edgey )
540  {
541  xerr = 31.0f * micronsToCm;
542  yerr = 90.0f * micronsToCm;
543  }
544  else
545  {
546  throw cms::Exception(" PixelCPETemplateReco::localError: Something wrong with pixel edge flag !!!");
547  }
548 
549  //cout << "xerr = " << xerr << endl;
550  //cout << "yerr = " << yerr << endl;
551  }
552  else
553  {
554  // &&& need a class const
555  //const float micronsToCm = 1.0e-4;
556 
557  xerr = theClusterParam.templSigmaX_ * micronsToCm;
558  yerr = theClusterParam.templSigmaY_ * micronsToCm;
559 
560  //cout << "xerr = " << xerr << endl;
561  //cout << "yerr = " << yerr << endl;
562 
563  // &&& should also check ierr (saved as class variable) and return
564  // &&& nonsense (another class static) if the template fit failed.
565  }
566 
567  if (theVerboseLevel > 9)
568  {
569  LogDebug("PixelCPETemplateReco") <<
570  " Sizex = " << theClusterParam.theCluster->sizeX() << " Sizey = " << theClusterParam.theCluster->sizeY() << " Edgex = " << edgex << " Edgey = " << edgey <<
571  " ErrX = " << xerr << " ErrY = " << yerr;
572  }
573 
574  } // else
575 
576  if ( !(xerr > 0.0f) )
577  throw cms::Exception("PixelCPETemplateReco::localError")
578  << "\nERROR: Negative pixel error xerr = " << xerr << "\n\n";
579 
580  if ( !(yerr > 0.0f) )
581  throw cms::Exception("PixelCPETemplateReco::localError")
582  << "\nERROR: Negative pixel error yerr = " << yerr << "\n\n";
583 
584  //cout << "Final errors set to: " << endl;
585  //cout << "xerr = " << xerr << endl;
586  //cout << "yerr = " << yerr << endl;
587  //cout << "Out of PixelCPETemplateREco..........................................................................." << endl;
588  //cout << endl;
589 
590  return LocalError(xerr*xerr, 0, yerr*yerr);
591 }
#define LogDebug(id)
bool isBarrel(GeomDetEnumerators::SubDetector m)
T sqrt(T t)
Definition: SSEVec.h:18
double f[11][100]
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
LocalPoint PixelCPETemplateReco::localPosition ( DetParam const &  theDetParam,
ClusterParam theClusterParam 
) const
privatevirtual

Implements PixelCPEBase.

Definition at line 116 of file PixelCPETemplateReco.cc.

References PixelCPEBase::DetParam::bz, gather_cfg::cout, PixelCPEBase::DetParam::detTemplateId, PixelCPEBase::DoLorentz_, GeomDet::geographicalId(), SiPixelTemplateDBObject::getTemplateID(), i, PixelCPETemplateReco::ClusterParamTemplate::ierr, GeomDetEnumerators::isEndcap(), RectangularPixelTopology::isItBigPixelInX(), RectangularPixelTopology::isItBigPixelInY(), GeomDetEnumerators::isTrackerPixel(), PixelCPEBase::LoadTemplatesFromDB_, Topology::localPosition(), PixelTopology::localX(), PixelTopology::localY(), LogDebug, PixelCPEBase::DetParam::lorentzShiftInCmX, PixelCPEBase::DetParam::lorentzShiftInCmY, SiPixelTemplate::lorxbias(), SiPixelTemplate::lorybias(), SiPixelTemplateReco::PixelTempReco2D(), SiPixelTemplateSplit::PixelTempSplit(), SiPixelTemplate2D::pushfile(), speed_, mathSSE::sqrt(), PixelCPEBase::templateDBobject_, PixelCPETemplateReco::ClusterParamTemplate::templProbQ_, PixelCPETemplateReco::ClusterParamTemplate::templProbX_, PixelCPETemplateReco::ClusterParamTemplate::templProbY_, PixelCPETemplateReco::ClusterParamTemplate::templQbin_, PixelCPETemplateReco::ClusterParamTemplate::templSigmaX_, PixelCPETemplateReco::ClusterParamTemplate::templSigmaY_, PixelCPETemplateReco::ClusterParamTemplate::templXrec_, PixelCPETemplateReco::ClusterParamTemplate::templYrec_, PixelCPEBase::DetParam::theDet, PixelCPEBase::DetParam::thePart, thePixelTemp_, PixelCPEBase::DetParam::theRecTopol, PixelCPEBase::DetParam::theTopol, create_public_lumi_plots::tmp_x, create_public_lumi_plots::tmp_y, unlikely, UseClusterSplitter_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

117 {
118 
119  ClusterParamTemplate & theClusterParam = static_cast<ClusterParamTemplate &>(theClusterParamBase);
120 
121  if(!GeomDetEnumerators::isTrackerPixel(theDetParam.thePart))
122  throw cms::Exception("PixelCPETemplateReco::localPosition :")
123  << "A non-pixel detector type in here?";
124  // barrel(false) or forward(true)
125  const bool fpix = GeomDetEnumerators::isEndcap(theDetParam.thePart);
126 
127  int ID = -9999;
128  if ( LoadTemplatesFromDB_ ) {
129  int ID0 = templateDBobject_->getTemplateID(theDetParam.theDet->geographicalId()); // just to comapre
130  ID = theDetParam.detTemplateId;
131  if(ID0!=ID) cout<<" different id"<< ID<<" "<<ID0<<endl;
132  } else { // from asci file
133  if ( !fpix ) ID = 40; // barrel
134  else ID = 41; // endcap
135  }
136  //cout << "PixelCPETemplateReco : ID = " << ID << endl;
137 
139 
140  // Make from cluster (a SiPixelCluster) a boost multi_array_2d called
141  // clust_array_2d.
142  boost::multi_array<float, 2> clust_array_2d(boost::extents[cluster_matrix_size_x][cluster_matrix_size_y]);
143 
144  // Preparing to retrieve ADC counts from the SiPixeltheClusterParam.theCluster-> In the cluster,
145  // we have the following:
146  // int minPixelRow(); // Minimum pixel index in the x direction (low edge).
147  // int maxPixelRow(); // Maximum pixel index in the x direction (top edge).
148  // int minPixelCol(); // Minimum pixel index in the y direction (left edge).
149  // int maxPixelCol(); // Maximum pixel index in the y direction (right edge).
150  // So the pixels from minPixelRow() will go into clust_array_2d[0][*],
151  // and the pixels from minPixelCol() will go into clust_array_2d[*][0].
152 
153  int row_offset = theClusterParam.theCluster->minPixelRow();
154  int col_offset = theClusterParam.theCluster->minPixelCol();
155 
156  // Store the coordinates of the center of the (0,0) pixel of the array that
157  // gets passed to PixelTempReco2D
158  // Will add these values to the output of PixelTempReco2D
159  float tmp_x = float(row_offset) + 0.5f;
160  float tmp_y = float(col_offset) + 0.5f;
161 
162  // Store these offsets (to be added later) in a LocalPoint after tranforming
163  // them from measurement units (pixel units) to local coordinates (cm)
164 
165  // ggiurgiu@jhu.edu 12/09/2010 : update call with trk angles needed for bow/kink corrections
166  LocalPoint lp;
167 
168  if ( theClusterParam.with_track_angle )
169  lp = theDetParam.theTopol->localPosition( MeasurementPoint(tmp_x, tmp_y), theClusterParam.loc_trk_pred );
170  else
171  {
172  edm::LogError("PixelCPETemplateReco")
173  << "@SUB = PixelCPETemplateReco::localPosition"
174  << "Should never be here. PixelCPETemplateReco should always be called with track angles. This is a bad error !!! ";
175 
176  lp = theDetParam.theTopol->localPosition( MeasurementPoint(tmp_x, tmp_y) );
177  }
178 
179  // Copy clust's pixels (calibrated in electrons) into clust_array_2d;
180  for (int i=0 ; i!=theClusterParam.theCluster->size(); ++i )
181  {
182  auto pix = theClusterParam.theCluster->pixel(i);
183  // *pixIter dereferences to Pixel struct, with public vars x, y, adc (all float)
184  // 02/13/2008 ggiurgiu@fnal.gov: type of x, y and adc has been changed to unsigned char, unsigned short, unsigned short
185  // in DataFormats/SiPixelCluster/interface/SiPixeltheClusterParam.theCluster->h so the type cast to int is redundant. Leave it there, it
186  // won't hurt.
187  int irow = int(pix.x) - row_offset; // &&& do we need +0.5 ???
188  int icol = int(pix.y) - col_offset; // &&& do we need +0.5 ???
189 
190  // Gavril : what do we do here if the row/column is larger than cluster_matrix_size_x/cluster_matrix_size_y = 7/21 ?
191  // Ignore them for the moment...
192  if ( irow<cluster_matrix_size_x && icol<cluster_matrix_size_y )
193  // 02/13/2008 ggiurgiu@fnal.gov typecast pixIter->adc to float
194  clust_array_2d[irow][icol] = (float)pix.adc;
195  }
196 
197  // Make and fill the bool arrays flagging double pixels
198  std::vector<bool> ydouble(cluster_matrix_size_y), xdouble(cluster_matrix_size_x);
199  // x directions (shorter), rows
200  for (int irow = 0; irow < cluster_matrix_size_x; ++irow)
201  {
202  xdouble[irow] = theDetParam.theRecTopol->isItBigPixelInX( irow+row_offset );
203  }
204 
205  // y directions (longer), columns
206  for (int icol = 0; icol < cluster_matrix_size_y; ++icol)
207  {
208  ydouble[icol] = theDetParam.theRecTopol->isItBigPixelInY( icol+col_offset );
209  }
210 
211  // Output:
212  float nonsense = -99999.9f; // nonsense init value
213  theClusterParam.templXrec_ = theClusterParam.templYrec_ = theClusterParam.templSigmaX_ = theClusterParam.templSigmaY_ = nonsense;
214  // If the template recontruction fails, we want to return 1.0 for now
215  theClusterParam.templProbY_ = theClusterParam.templProbX_ = theClusterParam.templProbQ_ = 1.0f;
216  theClusterParam.templQbin_ = 0;
217  // We have a boolean denoting whether the reco failed or not
218  theClusterParam.hasFilledProb_ = false;
219 
220  float templYrec1_ = nonsense;
221  float templXrec1_ = nonsense;
222  float templYrec2_ = nonsense;
223  float templXrec2_ = nonsense;
224 
225  // ******************************************************************
226  // Do it! Use cotalpha_ and cotbeta_ calculated in PixelCPEBase
227 
228 
229  float locBz = theDetParam.bz;
230 
231  theClusterParam.ierr =
232  PixelTempReco2D( ID, theClusterParam.cotalpha, theClusterParam.cotbeta,
233  locBz,
234  clust_array_2d, ydouble, xdouble,
235  templ,
236  theClusterParam.templYrec_, theClusterParam.templSigmaY_, theClusterParam.templProbY_,
237  theClusterParam.templXrec_, theClusterParam.templSigmaX_, theClusterParam.templProbX_,
238  theClusterParam.templQbin_,
239  speed_,
240  theClusterParam.templProbQ_
241  );
242 
243  // ******************************************************************
244 
245  // Check exit status
246  if unlikely( theClusterParam.ierr != 0 )
247  {
248  LogDebug("PixelCPETemplateReco::localPosition") <<
249  "reconstruction failed with error " << theClusterParam.ierr << "\n";
250 
251  // Gavril: what do we do in this case ? For now, just return the cluster center of gravity in microns
252  // In the x case, apply a rough Lorentz drift average correction
253  // To do: call PixelCPEGeneric whenever PixelTempReco2D fails
254  float lorentz_drift = -999.9;
255  if ( !fpix )
256  lorentz_drift = 60.0f; // in microns
257  else
258  lorentz_drift = 10.0f; // in microns
259  // ggiurgiu@jhu.edu, 21/09/2010 : trk angles needed to correct for bows/kinks
260  if ( theClusterParam.with_track_angle )
261  {
262  theClusterParam.templXrec_ = theDetParam.theTopol->localX( theClusterParam.theCluster->x(), theClusterParam.loc_trk_pred ) - lorentz_drift * micronsToCm; // rough Lorentz drift correction
263  theClusterParam.templYrec_ = theDetParam.theTopol->localY( theClusterParam.theCluster->y(), theClusterParam.loc_trk_pred );
264  }
265  else
266  {
267  edm::LogError("PixelCPETemplateReco")
268  << "@SUB = PixelCPETemplateReco::localPosition"
269  << "Should never be here. PixelCPETemplateReco should always be called with track angles. This is a bad error !!! ";
270 
271  theClusterParam.templXrec_ = theDetParam.theTopol->localX( theClusterParam.theCluster->x() ) - lorentz_drift * micronsToCm; // rough Lorentz drift correction
272  theClusterParam.templYrec_ = theDetParam.theTopol->localY( theClusterParam.theCluster->y() );
273  }
274  }
275  else if unlikely( UseClusterSplitter_ && theClusterParam.templQbin_ == 0 )
276  {
277  cout << " PixelCPETemplateReco : We should never be here !!!!!!!!!!!!!!!!!!!!!!" << endl;
278  cout << " (int)UseClusterSplitter_ = " << (int)UseClusterSplitter_ << endl;
279 
280  //ierr =
281  //PixelTempSplit( ID, fpix, cotalpha_, cotbeta_,
282  // clust_array_2d, ydouble, xdouble,
283  // templ,
284  // templYrec1_, templYrec2_, templSigmaY_, templProbY_,
285  // templXrec1_, templXrec2_, templSigmaX_, templProbX_,
286  // templQbin_ );
287 
288 
289  float dchisq;
290  float templProbQ_;
291  std::vector< SiPixelTemplateStore2D > thePixelTemp2D_;
292  SiPixelTemplate2D::pushfile(ID, thePixelTemp2D_);
293  SiPixelTemplate2D templ2D_(thePixelTemp2D_);
294 
295  theClusterParam.ierr =
296  SiPixelTemplateSplit::PixelTempSplit( ID, theClusterParam.cotalpha, theClusterParam.cotbeta,
297  clust_array_2d,
298  ydouble, xdouble,
299  templ,
300  templYrec1_, templYrec2_, theClusterParam.templSigmaY_, theClusterParam.templProbY_,
301  templXrec1_, templXrec2_, theClusterParam.templSigmaX_, theClusterParam.templProbX_,
302  theClusterParam.templQbin_,
303  templProbQ_,
304  true,
305  dchisq,
306  templ2D_ );
307 
308 
309  if ( theClusterParam.ierr != 0 )
310  {
311  LogDebug("PixelCPETemplateReco::localPosition") <<
312  "reconstruction failed with error " << theClusterParam.ierr << "\n";
313 
314  // Gavril: what do we do in this case ? For now, just return the cluster center of gravity in microns
315  // In the x case, apply a rough Lorentz drift average correction
316  // To do: call PixelCPEGeneric whenever PixelTempReco2D fails
317  float lorentz_drift = -999.9f;
318  if ( !fpix )
319  lorentz_drift = 60.0f; // in microns
320  else
321  lorentz_drift = 10.0f; // in microns
322 
323  // ggiurgiu@jhu.edu, 12/09/2010 : trk angles needed to correct for bows/kinks
324  if ( theClusterParam.with_track_angle )
325  {
326  theClusterParam.templXrec_ = theDetParam.theTopol->localX( theClusterParam.theCluster->x(),theClusterParam.loc_trk_pred ) - lorentz_drift * micronsToCm; // rough Lorentz drift correction
327  theClusterParam.templYrec_ = theDetParam.theTopol->localY( theClusterParam.theCluster->y(),theClusterParam.loc_trk_pred );
328  }
329  else
330  {
331  edm::LogError("PixelCPETemplateReco")
332  << "@SUB = PixelCPETemplateReco::localPosition"
333  << "Should never be here. PixelCPETemplateReco should always be called with track angles. This is a bad error !!! ";
334 
335  theClusterParam.templXrec_ = theDetParam.theTopol->localX( theClusterParam.theCluster->x() ) - lorentz_drift * micronsToCm; // very rough Lorentz drift correction
336  theClusterParam.templYrec_ = theDetParam.theTopol->localY( theClusterParam.theCluster->y() );
337 
338  }
339  }
340  else
341  {
342  // go from micrometer to centimeter
343  templXrec1_ *= micronsToCm;
344  templYrec1_ *= micronsToCm;
345  templXrec2_ *= micronsToCm;
346  templYrec2_ *= micronsToCm;
347 
348  // go back to the module coordinate system
349  templXrec1_ += lp.x();
350  templYrec1_ += lp.y();
351  templXrec2_ += lp.x();
352  templYrec2_ += lp.y();
353 
354  // calculate distance from each hit to the track and choose the
355  // hit closest to the track
356  float distance11 = sqrt( (templXrec1_ - theClusterParam.trk_lp_x)*(templXrec1_ - theClusterParam.trk_lp_x) +
357  (templYrec1_ - theClusterParam.trk_lp_y)*(templYrec1_ - theClusterParam.trk_lp_y) );
358 
359  float distance12 = sqrt( (templXrec1_ - theClusterParam.trk_lp_x)*(templXrec1_ - theClusterParam.trk_lp_x) +
360  (templYrec2_ - theClusterParam.trk_lp_y)*(templYrec2_ - theClusterParam.trk_lp_y) );
361 
362  float distance21 = sqrt( (templXrec2_ - theClusterParam.trk_lp_x)*(templXrec2_ - theClusterParam.trk_lp_x) +
363  (templYrec1_ - theClusterParam.trk_lp_y)*(templYrec1_ - theClusterParam.trk_lp_y) );
364 
365  float distance22 = sqrt( (templXrec2_ - theClusterParam.trk_lp_x)*(templXrec2_ - theClusterParam.trk_lp_x) +
366  (templYrec2_ - theClusterParam.trk_lp_y)*(templYrec2_ - theClusterParam.trk_lp_y) );
367 
368  float min_templXrec_ = -999.9;
369  float min_templYrec_ = -999.9;
370  float distance_min = 9999999999.9;
371  if ( distance11 < distance_min )
372  {
373  distance_min = distance11;
374  min_templXrec_ = templXrec1_;
375  min_templYrec_ = templYrec1_;
376  }
377  if ( distance12 < distance_min )
378  {
379  distance_min = distance12;
380  min_templXrec_ = templXrec1_;
381  min_templYrec_ = templYrec2_;
382  }
383  if ( distance21 < distance_min )
384  {
385  distance_min = distance21;
386  min_templXrec_ = templXrec2_;
387  min_templYrec_ = templYrec1_;
388  }
389  if ( distance22 < distance_min )
390  {
391  distance_min = distance22;
392  min_templXrec_ = templXrec2_;
393  min_templYrec_ = templYrec2_;
394  }
395 
396  theClusterParam.templXrec_ = min_templXrec_;
397  theClusterParam.templYrec_ = min_templYrec_;
398  }
399  } // else if ( UseClusterSplitter_ && templQbin_ == 0 )
400 
401  else // apparenly this is he good one!
402  {
403  // go from micrometer to centimeter
404  theClusterParam.templXrec_ *= micronsToCm;
405  theClusterParam.templYrec_ *= micronsToCm;
406 
407  // go back to the module coordinate system
408  theClusterParam.templXrec_ += lp.x();
409  theClusterParam.templYrec_ += lp.y();
410 
411  // Compute the Alignment Group Corrections [template ID should already be selected from call to reco procedure]
412  if ( DoLorentz_ ) {
413  // Do only if the lotentzshift has meaningfull numbers
414  if( theDetParam.lorentzShiftInCmX!= 0.0 || theDetParam.lorentzShiftInCmY!= 0.0 ) {
415  // the LA width/shift returned by templates use (+)
416  // the LA width/shift produced by PixelCPEBase for positive LA is (-)
417  // correct this by iserting (-)
418  //float temp1 = -micronsToCm*templ.lorxwidth(); // old
419  //float temp2 = -micronsToCm*templ.lorywidth(); // does not incl 1/2
420  float templateLorbiasCmX = -micronsToCm*templ.lorxbias(); // new
421  float templateLorbiasCmY = -micronsToCm*templ.lorybias(); //incl. 1/2
422  // now, correctly, we can use the difference of shifts
423  //theClusterParam.templXrec_ += 0.5*(theDetParam.lorentzShiftInCmX - templateLorbiasCmX);
424  //theClusterParam.templYrec_ += 0.5*(theDetParam.lorentzShiftInCmY - templateLorbiasCmY);
425  theClusterParam.templXrec_ += (0.5*(theDetParam.lorentzShiftInCmX) - templateLorbiasCmX);
426  theClusterParam.templYrec_ += (0.5*(theDetParam.lorentzShiftInCmY) - templateLorbiasCmY);
427  //cout << "Templates: la lorentz offset = "
428  // <<(0.5*(theDetParam.lorentzShiftInCmX)-templateLorbiasCmX)
429  // <<" "<<templateLorbiasCmX<<" "<<templateLorbiasCmY
430  // <<" "<<temp1<<" "<<temp2
431  // <<" "<<theDetParam.lorentzShiftInCmX
432  // <<" "<<theDetParam.lorentzShiftInCmY
433  // << endl; //dk
434  } //else {cout<<" LA is 0, disable offset corrections "<<endl;} //dk
435  } //else {cout<<" Do not do LA offset correction "<<endl;} //dk
436 
437  }
438 
439  // Save probabilities and qBin in the quantities given to us by the base class
440  // (for which there are also inline getters). &&& templProbX_ etc. should be retired...
441  theClusterParam.probabilityX_ = theClusterParam.templProbX_;
442  theClusterParam.probabilityY_ = theClusterParam.templProbY_;
443  theClusterParam.probabilityQ_ = theClusterParam.templProbQ_;
444  theClusterParam.qBin_ = theClusterParam.templQbin_;
445 
446  if ( theClusterParam.ierr == 0 ) // always true here
447  theClusterParam.hasFilledProb_ = true;
448 
449  return LocalPoint( theClusterParam.templXrec_, theClusterParam.templYrec_ );
450 
451 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
short getTemplateID(const uint32_t &detid) const
static bool pushfile(int filenum, std::vector< SiPixelTemplateStore2D > &thePixelTemp_)
uint32_t ID
Definition: Definitions.h:26
T y() const
Definition: PV3DBase.h:63
std::vector< SiPixelTemplateStore > thePixelTemp_
#define unlikely(x)
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:255
T sqrt(T t)
Definition: SSEVec.h:18
int PixelTempSplit(int id, float cotalpha, float cotbeta, array_2d &cluster, std::vector< bool > &ydouble, std::vector< bool > &xdouble, SiPixelTemplate &templ, float &yrec1, float &yrec2, float &sigmay, float &prob2y, float &xrec1, float &xrec2, float &sigmax, float &prob2x, int &q2bin, float &prob2Q, bool resolve, int speed, float &dchisq, bool deadpix, std::vector< std::pair< int, int > > &zeropix, SiPixelTemplate2D &templ2D)
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:251
int PixelTempReco2D(int id, float cotalpha, float cotbeta, float locBz, array_2d &cluster, std::vector< bool > &ydouble, std::vector< bool > &xdouble, 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)
bool isEndcap(GeomDetEnumerators::SubDetector m)
bool isTrackerPixel(const GeomDetEnumerators::SubDetector m)
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
tuple cout
Definition: gather_cfg.py:145
T x() const
Definition: PV3DBase.h:62

Member Data Documentation

int PixelCPETemplateReco::speed_
private

Definition at line 78 of file PixelCPETemplateReco.h.

Referenced by localPosition(), and PixelCPETemplateReco().

std::vector< SiPixelTemplateStore > PixelCPETemplateReco::thePixelTemp_
private

Definition at line 76 of file PixelCPETemplateReco.h.

Referenced by localPosition(), and PixelCPETemplateReco().

bool PixelCPETemplateReco::UseClusterSplitter_
private

Definition at line 80 of file PixelCPETemplateReco.h.

Referenced by localPosition(), and PixelCPETemplateReco().