CMS 3D CMS Logo

PixelCPEBase.cc
Go to the documentation of this file.
1 // Move geomCorrection to the concrete class. d.k. 06/06.
2 // Change drift direction. d.k. 06/06
3 // G. Giurgiu (ggiurgiu@pha.jhu.edu), 12/01/06, implemented the function:
4 // computeAnglesFromDetPosition(const SiPixelCluster & cl,
5 // change to use Lorentz angle from DB Lotte Wilke, Jan. 31st, 2008
6 // Change to use Generic error & Template calibration from DB - D.Fehling 11/08
7 
8 
12 
14 
15 #define CORRECT_FOR_BIG_PIXELS
16 
17 // MessageLogger
19 
20 // Magnetic field
22 
23 
24 #include <iostream>
25 
26 using namespace std;
27 
28 //-----------------------------------------------------------------------------
29 // A constructor run for generic and templates
30 //
31 //-----------------------------------------------------------------------------
33  const MagneticField *mag,
34  const TrackerGeometry& geom,
35  const TrackerTopology& ttopo,
36  const SiPixelLorentzAngle * lorentzAngle,
37  const SiPixelGenErrorDBObject * genErrorDBObject,
38  const SiPixelTemplateDBObject * templateDBobject,
39  const SiPixelLorentzAngle * lorentzAngleWidth,
40  int flag)
41 // : useLAAlignmentOffsets_(false), useLAOffsetFromConfig_(false),
42 : useLAOffsetFromConfig_(false),
43 useLAWidthFromConfig_(false), useLAWidthFromDB_(false), theFlag_(flag),
44 magfield_(mag), geom_(geom), ttopo_(ttopo)
45 {
46 
47 #ifdef EDM_ML_DEBUG
48  nRecHitsTotal_=0;
49  nRecHitsUsedEdge_=0,
50 #endif
51 
52  //--- Lorentz angle tangent per Tesla
53  lorentzAngle_ = lorentzAngle;
54  lorentzAngleWidth_ = lorentzAngleWidth;
55 
56  //-- GenError Calibration Object (different from SiPixelCPEGenericErrorParm) from DB
57  genErrorDBObject_ = genErrorDBObject;
58  //cout<<" new errors "<<genErrorDBObject<<" "<<genErrorDBObject_<<endl;
59 
60  //-- Template Calibration Object from DB
61  if(theFlag_!=0) templateDBobject_ = templateDBobject; // flag to check if it is generic or templates
62 
63  // Configurables
64  // For both templates & generic
65 
66  // Read templates and/or generic errors from DB
67  LoadTemplatesFromDB_ = conf.getParameter<bool>("LoadTemplatesFromDB");
68  //cout<<" use generros/templaets "<<LoadTemplatesFromDB_<<endl;
69 
70  //--- Algorithm's verbosity
72  conf.getUntrackedParameter<int>("VerboseLevel",0);
73 
74  //-- Switch on/off E.B
75  alpha2Order = conf.getParameter<bool>("Alpha2Order");
76 
77  //--- A flag that could be used to change the behavior of
78  //--- clusterProbability() in TSiPixelRecHit (the *transient* one).
79  //--- The problem is that the transient hits are made after the CPE runs
80  //--- and they don't get the access to the PSet, so we pass it via the
81  //--- CPE itself...
82  //
84  = (unsigned int) conf.getParameter<int>("ClusterProbComputationFlag");
85 
86  // This LA related parameters are only relevant for the Generic algo
87  // They still have to be used in Base since the LA computation is in Base
88 
89  // Use LA-width from DB.
90  // If both (this and from config) are false LA-width is calcuated from LA-offset
91  useLAWidthFromDB_ = conf.existsAs<bool>("useLAWidthFromDB")?
92  conf.getParameter<bool>("useLAWidthFromDB"):false;
93 
94  // Use Alignment LA-offset in generic
95  //useLAAlignmentOffsets_ = conf.existsAs<bool>("useLAAlignmentOffsets")?
96  //conf.getParameter<bool>("useLAAlignmentOffsets"):false;
97 
98  // Used only for testing
99  lAOffset_ = conf.existsAs<double>("lAOffset")? // fixed LA value
100  conf.getParameter<double>("lAOffset"):0.0;
101  lAWidthBPix_ = conf.existsAs<double>("lAWidthBPix")? // fixed LA width
102  conf.getParameter<double>("lAWidthBPix"):0.0;
103  lAWidthFPix_ = conf.existsAs<double>("lAWidthFPix")? // fixed LA width
104  conf.getParameter<double>("lAWidthFPix"):0.0;
105 
106  // Use LA-offset from config, for testing only
107  if(lAOffset_>0.0) useLAOffsetFromConfig_ = true;
108  // Use LA-width from config, split into fpix & bpix, for testing only
109  if(lAWidthBPix_>0.0 || lAWidthFPix_>0.0) useLAWidthFromConfig_ = true;
110 
111 
112  // For Templates only
113  // Compute the Lorentz shifts for this detector element for templates (from Alignment)
114  DoLorentz_ = conf.existsAs<bool>("DoLorentz")?conf.getParameter<bool>("DoLorentz"):false;
115 
116  LogDebug("PixelCPEBase") <<" LA constants - "
117  <<lAOffset_<<" "<<lAWidthBPix_<<" "<<lAWidthFPix_<<endl; //dk
118 
119  fillDetParams();
120 
121  //cout<<" LA "<<lAOffset_<<" "<<lAWidthBPix_<<" "<<lAWidthFPix_<<endl; //dk
122 }
123 
124 //-----------------------------------------------------------------------------
125 // Fill all variables which are constant for an event (geometry)
126 //-----------------------------------------------------------------------------
128 {
129  //cout<<" in fillDetParams "<<theFlag_<<endl;
130 
131  auto const & dus = geom_.detUnits();
132  unsigned m_detectors = dus.size();
133  for(unsigned int i=1;i<7;++i) {
134  LogDebug("LookingForFirstStrip") << "Subdetector " << i
135  << " GeomDetEnumerator " << GeomDetEnumerators::tkDetEnum[i]
136  << " offset " << geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i])
137  << " is it strip? " << (geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) != dus.size() ?
138  dus[geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i])]->type().isTrackerStrip() : false);
139  if(geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i]) != dus.size() &&
140  dus[geom_.offsetDU(GeomDetEnumerators::tkDetEnum[i])]->type().isTrackerStrip()) {
142  }
143  }
144  LogDebug("LookingForFirstStrip") << " Chosen offset: " << m_detectors;
145 
146 
147  m_DetParams.resize(m_detectors);
148  //cout<<"caching "<<m_detectors<<" pixel detectors"<<endl;
149  for (unsigned i=0; i!=m_detectors;++i) {
150  auto & p=m_DetParams[i];
151  p.theDet = dynamic_cast<const PixelGeomDetUnit*>(dus[i]);
152  assert(p.theDet);
153  assert(p.theDet->index()==int(i));
154 
155  p.theOrigin = p.theDet->surface().toLocal(GlobalPoint(0,0,0));
156 
157  //--- p.theDet->type() returns a GeomDetType, which implements subDetector()
158  p.thePart = p.theDet->type().subDetector();
159 
160  //cout<<" in PixelCPEBase - in det "<<thePart<<endl; //dk
161 
162  //--- The location in of this DetUnit in a cyllindrical coord system (R,Z)
163  //--- The call goes via BoundSurface, returned by p.theDet->surface(), but
164  //--- position() is implemented in GloballyPositioned<> template
165  //--- ( BoundSurface : Surface : GloballyPositioned<float> )
166  //p.theDetR = p.theDet->surface().position().perp(); //Not used, AH
167  //p.theDetZ = p.theDet->surface().position().z(); //Not used, AH
168  //--- Define parameters for chargewidth calculation
169 
170  //--- bounds() is implemented in BoundSurface itself.
171  p.theThickness = p.theDet->surface().bounds().thickness();
172 
173  // Cache the det id for templates and generic erros
174 
175  if(theFlag_==0) { // for generic
176  if(LoadTemplatesFromDB_ ) // do only if genError requested
177  p.detTemplateId = genErrorDBObject_->getGenErrorID(p.theDet->geographicalId().rawId());
178  } else { // for templates
179  p.detTemplateId = templateDBobject_->getTemplateID(p.theDet->geographicalId());
180  }
181 
182  // just for testing
183  //int i1 = 0;
184  //if(theFlag_==0) i1 = genErrorDBObject_->getGenErrorID(p.theDet->geographicalId().rawId());
185  //int i2= templateDBobject_->getTemplateID(p.theDet->geographicalId().rawId());
186  //int i3= templateDBobject_->getTemplateID(p.theDet->geographicalId());
187  //if(i2!=i3) cout<<i2<<" != "<<i3<<endl;
188  //cout<<i<<" "<<p.detTemplateId<<" "<<i1<<" "<<i2<<" "<<i3<<endl;
189 
190  auto topol = &(p.theDet->specificTopology());
191  p.theTopol=topol;
192  auto const proxyT = dynamic_cast<const ProxyPixelTopology*>(p.theTopol);
193  if (proxyT) p.theRecTopol = dynamic_cast<const RectangularPixelTopology*>(&(proxyT->specificTopology()));
194  else p.theRecTopol = dynamic_cast<const RectangularPixelTopology*>(p.theTopol);
195  assert(p.theRecTopol);
196 
197  //---- The geometrical description of one module/plaquette
198  //p.theNumOfRow = p.theRecTopol->nrows(); // rows in x //Not used, AH
199  //p.theNumOfCol = p.theRecTopol->ncolumns(); // cols in y //Not used, AH
200  std::pair<float,float> pitchxy = p.theRecTopol->pitch();
201  p.thePitchX = pitchxy.first; // pitch along x
202  p.thePitchY = pitchxy.second; // pitch along y
203 
204  //p.theSign = isFlipped(&p) ? -1 : 1; //Not used, AH
205 
206  LocalVector Bfield = p.theDet->surface().toLocal(magfield_->inTesla(p.theDet->surface().position()));
207  p.bz = Bfield.z();
208  p.bx = Bfield.x();
209 
210 
211  // Compute the Lorentz shifts for this detector element
212  if ( (theFlag_==0) || DoLorentz_ ) { // do always for generic and if(DOLorentz) for templates
213  p.driftDirection = driftDirection(p, Bfield );
215  }
216 
217 
218  LogDebug("PixelCPEBase") << "***** PIXEL LAYOUT *****"
219  << " thePart = " << p.thePart
220  << " theThickness = " << p.theThickness
221  << " thePitchX = " << p.thePitchX
222  << " thePitchY = " << p.thePitchY;
223  // << " theLShiftX = " << p.theLShiftX;
224 
225 
226  }
227 }
228 
229 //-----------------------------------------------------------------------------
230 // One function to cache the variables common for one DetUnit.
231 //-----------------------------------------------------------------------------
232 void
233 PixelCPEBase::setTheClu( DetParam const & theDetParam, ClusterParam & theClusterParam ) const
234 {
235 
236  //--- Geometric Quality Information
237  int minInX,minInY,maxInX,maxInY=0;
238  minInX = theClusterParam.theCluster->minPixelRow();
239  minInY = theClusterParam.theCluster->minPixelCol();
240  maxInX = theClusterParam.theCluster->maxPixelRow();
241  maxInY = theClusterParam.theCluster->maxPixelCol();
242 
243  theClusterParam.isOnEdge_ = theDetParam.theRecTopol->isItEdgePixelInX(minInX) | theDetParam.theRecTopol->isItEdgePixelInX(maxInX) |
244  theDetParam.theRecTopol->isItEdgePixelInY(minInY) | theDetParam.theRecTopol->isItEdgePixelInY(maxInY) ;
245 
246  // FOR NOW UNUSED. KEEP IT IN CASE WE WANT TO USE IT IN THE FUTURE
247  // Bad Pixels have their charge set to 0 in the clusterizer
248  //hasBadPixels_ = false;
249  //for(unsigned int i=0; i<theClusterParam.theCluster->pixelADC().size(); ++i) {
250  //if(theClusterParam.theCluster->pixelADC()[i] == 0) { hasBadPixels_ = true; break;}
251  //}
252 
253  theClusterParam.spansTwoROCs_ = theDetParam.theRecTopol->containsBigPixelInX(minInX,maxInX) |
254  theDetParam.theRecTopol->containsBigPixelInY(minInY,maxInY);
255 
256 }
257 
258 
259 //-----------------------------------------------------------------------------
260 // Compute alpha_ and beta_ from the LocalTrajectoryParameters.
261 // Note: should become const after both localParameters() become const.
262 //-----------------------------------------------------------------------------
263 void PixelCPEBase::
264 computeAnglesFromTrajectory( DetParam const & theDetParam, ClusterParam & theClusterParam,
265  const LocalTrajectoryParameters & ltp) const
266 {
267  //cout<<" in PixelCPEBase:computeAnglesFromTrajectory - "<<endl; //dk
268  /*
269  //theClusterParam.loc_traj_param = ltp;
270  LocalVector localDir = ltp.momentum();
271  float locx = localDir.x();
272  float locy = localDir.y();
273  float locz = localDir.z();
274  // Danek's definition
275  alpha_ = acos(locx/sqrt(locx*locx+locz*locz));
276  if ( isFlipped() ) // &&& check for FPIX !!!
277  alpha_ = PI - alpha_ ;
278  beta_ = acos(locy/sqrt(locy*locy+locz*locz));
279  */
280 
281 
282  theClusterParam.cotalpha = ltp.dxdz();
283  theClusterParam.cotbeta = ltp.dydz();
284 
285 
286  LocalPoint trk_lp = ltp.position();
287  theClusterParam.trk_lp_x = trk_lp.x();
288  theClusterParam.trk_lp_y = trk_lp.y();
289 
290  theClusterParam.with_track_angle = true;
291 
292  // ggiurgiu@jhu.edu 12/09/2010 : needed to correct for bows/kinks
293  theClusterParam.loc_trk_pred =
294  Topology::LocalTrackPred(theClusterParam.trk_lp_x, theClusterParam.trk_lp_y,
295  theClusterParam.cotalpha,theClusterParam.cotbeta);
296 
297 }
298 
299 //-----------------------------------------------------------------------------
300 // Estimate theAlpha for barrel, based on the det position.
301 // &&& Needs to be consolidated from the above.
302 //-----------------------------------------------------------------------------
303 //float
304 //PixelCPEBase::estimatedAlphaForBarrel(float centerx) const
305 //{
306 // float tanalpha = theSign * (centerx-theOffsetX) * thePitchX / theDetR;
307 // return PI/2.0 - atan(tanalpha);
308 //}
309 
310 
311 
312 
313 //-----------------------------------------------------------------------------
314 // Compute alpha_ and beta_ from the position of this DetUnit.
315 // &&& DOES NOT WORK FOR NOW. d.k. 6/06
316 // The angles from dets are calculated internaly in the PixelCPEInitial class
317 //-----------------------------------------------------------------------------
318 // G. Giurgiu, 12/01/06 : implement the function
319 void PixelCPEBase::
320 computeAnglesFromDetPosition(DetParam const & theDetParam, ClusterParam & theClusterParam ) const
321 {
322 
323 
324  /*
325  // get cluster center of gravity (of charge)
326  float xcenter = cl.x();
327  float ycenter = cl.y();
328 
329  // get the cluster position in local coordinates (cm)
330 
331  // ggiurgiu@jhu.edu 12/09/2010 : This function is called without track info, therefore there are no track
332  // angles to provide here. Call the default localPosition (without track info)
333  LocalPoint lp = theTopol->localPosition( MeasurementPoint(xcenter, ycenter) );
334 
335 
336  // get the cluster position in global coordinates (cm)
337  GlobalPoint gp = theDet->surface().toGlobal( lp );
338  float gp_mod = sqrt( gp.x()*gp.x() + gp.y()*gp.y() + gp.z()*gp.z() );
339 
340  // normalize
341  float gpx = gp.x()/gp_mod;
342  float gpy = gp.y()/gp_mod;
343  float gpz = gp.z()/gp_mod;
344 
345  // make a global vector out of the global point; this vector will point from the
346  // origin of the detector to the cluster
347  GlobalVector gv(gpx, gpy, gpz);
348 
349  // make local unit vector along local X axis
350  const Local3DVector lvx(1.0, 0.0, 0.0);
351 
352  // get the unit X vector in global coordinates/
353  GlobalVector gvx = theDet->surface().toGlobal( lvx );
354 
355  // make local unit vector along local Y axis
356  const Local3DVector lvy(0.0, 1.0, 0.0);
357 
358  // get the unit Y vector in global coordinates
359  GlobalVector gvy = theDet->surface().toGlobal( lvy );
360 
361  // make local unit vector along local Z axis
362  const Local3DVector lvz(0.0, 0.0, 1.0);
363 
364  // get the unit Z vector in global coordinates
365  GlobalVector gvz = theDet->surface().toGlobal( lvz );
366 
367  // calculate the components of gv (the unit vector pointing to the cluster)
368  // in the local coordinate system given by the basis {gvx, gvy, gvz}
369  // note that both gv and the basis {gvx, gvy, gvz} are given in global coordinates
370  float gv_dot_gvx = gv.x()*gvx.x() + gv.y()*gvx.y() + gv.z()*gvx.z();
371  float gv_dot_gvy = gv.x()*gvy.x() + gv.y()*gvy.y() + gv.z()*gvy.z();
372  float gv_dot_gvz = gv.x()*gvz.x() + gv.y()*gvz.y() + gv.z()*gvz.z();
373  */
374 
375  // all the above is equivalent to
376  LocalPoint lp = theDetParam.theTopol->localPosition( MeasurementPoint(theClusterParam.theCluster->x(), theClusterParam.theCluster->y()) );
377  auto gvx = lp.x()-theDetParam.theOrigin.x();
378  auto gvy = lp.y()-theDetParam.theOrigin.y();
379  auto gvz = -1.f/theDetParam.theOrigin.z();
380  // normalization not required as only ratio used...
381 
382 
383  //theClusterParam.zneg = (gvz < 0); // Not used, AH
384 
385  // calculate angles
386  theClusterParam.cotalpha = gvx*gvz;
387  theClusterParam.cotbeta = gvy*gvz;
388 
389  theClusterParam.with_track_angle = false;
390 
391 
392  /*
393  // used only in dberror param...
394  auto alpha = HALF_PI - std::atan(cotalpha_);
395  auto beta = HALF_PI - std::atan(cotbeta_);
396  if (zneg) { beta -=PI; alpha -=PI;}
397 
398  auto alpha_ = atan2( gv_dot_gvz, gv_dot_gvx );
399  auto beta_ = atan2( gv_dot_gvz, gv_dot_gvy );
400 
401  assert(std::abs(std::round(alpha*10000.f)-std::round(alpha_*10000.f))<2);
402  assert(std::abs(std::round(beta*10000.f)-std::round(beta_*10000.f))<2);
403  */
404 
405 }
406 
407 
408 
409 //-----------------------------------------------------------------------------
410 // The isFlipped() is a silly way to determine which detectors are inverted.
411 // In the barrel for every 2nd ladder the E field direction is in the
412 // global r direction (points outside from the z axis), every other
413 // ladder has the E field inside. Something similar is in the
414 // forward disks (2 sides of the blade). This has to be recognised
415 // because the charge sharing effect is different.
416 //
417 // The isFliped does it by looking and the relation of the local (z always
418 // in the E direction) to global coordinates. There is probably a much
419 // better way.
420 //-----------------------------------------------------------------------------
421 bool PixelCPEBase::isFlipped(DetParam const & theDetParam) const
422 {
423  // Check the relative position of the local +/- z in global coordinates.
424  float tmp1 = theDetParam.theDet->surface().toGlobal(Local3DPoint(0.,0.,0.)).perp2();
425  float tmp2 = theDetParam.theDet->surface().toGlobal(Local3DPoint(0.,0.,1.)).perp2();
426  //cout << " 1: " << tmp1 << " 2: " << tmp2 << endl;
427  if ( tmp2<tmp1 ) return true;
428  else return false;
429 }
430 //------------------------------------------------------------------------
432  auto i = det.index();
433  //cout << "get parameters of detector " << i << endl;
434  assert(i<int(m_DetParams.size()));
435  //if (i>=int(m_DetParams.size())) m_DetParams.resize(i+1); // should never happen!
436  const DetParam & p = m_DetParams[i];
437  return p;
438 }
439 
440 //-----------------------------------------------------------------------------
441 // Drift direction.
442 // Works OK for barrel and forward.
443 // The formulas used for dir_x,y,z have to be exactly the same as the ones
444 // used in the digitizer (SiPixelDigitizerAlgorithm.cc).
445 //
446 //-----------------------------------------------------------------------------
448 PixelCPEBase::driftDirection(DetParam & theDetParam, GlobalVector bfield ) const {
449 
450  Frame detFrame(theDetParam.theDet->surface().position(), theDetParam.theDet->surface().rotation());
451  LocalVector Bfield = detFrame.toLocal(bfield);
452  return driftDirection(theDetParam,Bfield);
453 
454 }
455 
457 PixelCPEBase::driftDirection(DetParam & theDetParam, LocalVector Bfield ) const {
458  const bool LocalPrint = false;
459 
460  // Use LA from DB or from config
461  float langle = 0.;
462  if( !useLAOffsetFromConfig_ ) { // get it from DB
463  if(lorentzAngle_ != nullptr) { // a real LA object
464  langle = lorentzAngle_->getLorentzAngle(theDetParam.theDet->geographicalId().rawId());
465  //cout<<" la "<<langle<<" "<< theDetParam.theDet->geographicalId().rawId() <<endl;
466  } else { // no LA, unused
467  //cout<<" LA object is NULL, assume LA = 0"<<endl; //dk
468  langle = 0; // set to a fake value
469  }
470  if(LocalPrint) cout<<" Will use LA Offset from DB "<<langle<<endl;
471  } else { // from config file
472  langle = lAOffset_;
473  if(LocalPrint) cout<<" Will use LA Offset from config "<<langle<<endl;
474  }
475 
476  // Now do the LA width stuff
477  theDetParam.widthLAFractionX = 1.; // predefine to 1 (default) if things fail
478  theDetParam.widthLAFractionY = 1.;
479 
480  // Compute the charge width, generic only
481  if(theFlag_==0) {
482 
483  if(useLAWidthFromDB_ && (lorentzAngleWidth_ != nullptr) ) {
484  // take it from a seperate, special LA DB object (forWidth)
485 
486  auto langleWidth = lorentzAngleWidth_->getLorentzAngle(theDetParam.theDet->geographicalId().rawId());
487  if(langleWidth!=0.0) theDetParam.widthLAFractionX = std::abs(langleWidth/langle);
488  // leave the widthLAFractionY=1.
489  //cout<<" LAWidth lorentz width "<<theDetParam.widthLAFractionX<<" "<<theDetParam.widthLAFractionY<<endl;
490 
491  } else if(useLAWidthFromConfig_) { // get from config
492 
493  double lAWidth=0;
494  if( GeomDetEnumerators::isTrackerPixel(theDetParam.thePart) && GeomDetEnumerators::isBarrel(theDetParam.thePart) ) lAWidth = lAWidthBPix_; // barrel
495  else lAWidth = lAWidthFPix_;
496 
497  if(langle!=0.0) theDetParam.widthLAFractionX = std::abs(lAWidth/langle);
498  // fix the FractionY at 1
499 
500  //cout<<" Lorentz width from config "<<theDetParam.widthLAFractionX<<" "<<theDetParam.widthLAFractionY<<endl;
501 
502  } else { // get if from the offset LA (old method used until 2013)
503  // do nothing
504  //cout<<" Old default LA width method "<<theDetParam.widthLAFractionX<<" "<<theDetParam.widthLAFractionY<<endl;
505 
506  }
507 
508  //cout<<" Final LA fraction "<<theDetParam.widthLAFractionX<<" "<<theDetParam.widthLAFractionY<<endl;
509 
510  } // if flag
511 
512 
513  if(LocalPrint) cout<<" in PixelCPEBase:driftDirection - "<<langle<<" "<<Bfield<<endl; //dk
514 
515  float alpha2 = alpha2Order ? langle*langle : 0; //
516 
517  // **********************************************************************
518  // Our convention is the following:
519  // +x is defined by the direction of the Lorentz drift!
520  // +z is defined by the direction of E field (so electrons always go into -z!)
521  // +y is defined by +x and +z, and it turns out to be always opposite to the +B field.
522  // **********************************************************************
523 
524  float dir_x = -( langle * Bfield.y() + alpha2* Bfield.z()* Bfield.x() );
525  float dir_y = ( langle * Bfield.x() - alpha2* Bfield.z()* Bfield.y() );
526  float dir_z = -( 1.f + alpha2* Bfield.z()*Bfield.z() );
527  auto scale = 1.f/std::abs( dir_z ); // same as 1 + alpha2*Bfield.z()*Bfield.z()
528  LocalVector dd(dir_x*scale, dir_y*scale, -1.f ); // last is -1 !
529 
530  LogDebug("PixelCPEBase") << " The drift direction in local coordinate is " << dd ;
531 
532  return dd;
533 }
534 
535 //-----------------------------------------------------------------------------
536 // One-shot computation of the driftDirection and both lorentz shifts
537 //-----------------------------------------------------------------------------
538 void
540 
541  //cout<<" in PixelCPEBase:computeLorentzShifts - "<<driftDirection_<<endl; //dk
542 
543  // Max shift (at the other side of the sensor) in cm
544  theDetParam.lorentzShiftInCmX = theDetParam.driftDirection.x()/theDetParam.driftDirection.z() * theDetParam.theThickness; //
545  theDetParam.lorentzShiftInCmY = theDetParam.driftDirection.y()/theDetParam.driftDirection.z() * theDetParam.theThickness; //
546 
547  //cout<<" in PixelCPEBase:computeLorentzShifts - "
548  //<<lorentzShiftInCmX_<<" "
549  //<<lorentzShiftInCmY_<<" "
550  //<<endl; //dk
551 
552  LogDebug("PixelCPEBase") << " The drift direction in local coordinate is "
553  << theDetParam.driftDirection ;
554 }
555 
556 //-----------------------------------------------------------------------------
561 //-----------------------------------------------------------------------------
564 {
566  if (theClusterParam.hasFilledProb_) {
567  float probabilityXY=0;
568  if ( theClusterParam.probabilityX_ !=0 && theClusterParam.probabilityY_ !=0 )
569  probabilityXY = theClusterParam.probabilityX_ * theClusterParam.probabilityY_ * (1.f - std::log(theClusterParam.probabilityX_ * theClusterParam.probabilityY_) ) ;
571  qualWord );
572 
574  qualWord );
575  }
577  qualWord );
578 
580  qualWord );
581 
583  qualWord );
584 
586  qualWord );
587 
589  qualWord );
590 
591  return qualWord;
592 }
#define LogDebug(id)
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
void setProbabilityXY(float prob, QualWordType &qualWord) const
T getUntrackedParameter(std::string const &, T const &) const
static const Packing thePacking
int minPixelCol() const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
Local3DPoint theOrigin
Definition: PixelCPEBase.h:66
const SiPixelLorentzAngle * lorentzAngleWidth_
Definition: PixelCPEBase.h:244
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
LocalVector driftDirection(DetParam &theDetParam, GlobalVector bfield) const
short getTemplateID(const uint32_t &detid) const
LocalPoint position() const
Local x and y position coordinates.
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:87
bool isBarrel(GeomDetEnumerators::SubDetector m)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
bool isFlipped(DetParam const &theDetParam) const
DetParams m_DetParams
Definition: PixelCPEBase.h:280
T perp2() const
Definition: PV3DBase.h:71
void fillDetParams()
short getGenErrorID(const uint32_t &detid) const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:60
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:65
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:63
bool useLAWidthFromDB_
Definition: PixelCPEBase.h:233
float lAWidthBPix_
Definition: PixelCPEBase.h:228
int maxPixelRow() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
bool useLAWidthFromConfig_
Definition: PixelCPEBase.h:232
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
bool isItEdgePixelInX(int ixbin) const override
const PixelTopology * theTopol
Definition: PixelCPEBase.h:62
int minPixelRow() const
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:253
void computeAnglesFromTrajectory(DetParam const &theDetParam, ClusterParam &theClusterParam, const LocalTrajectoryParameters &ltp) const
unsigned int offsetDU(SubDetector sid) const
T z() const
Definition: PV3DBase.h:64
void computeLorentzShifts(DetParam &) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:249
void setHasBadPixels(bool flag, QualWordType &qualWord) const
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
double f[11][100]
const MagneticField * magfield_
Definition: PixelCPEBase.h:239
int index() const
Definition: GeomDet.h:99
void computeAnglesFromDetPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const
void setProbabilityQ(float prob, QualWordType &qualWord) const
SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam &theClusterParam) const
SubDetector tkDetEnum[8]
void setSpansTwoROCs(bool flag, QualWordType &qualWord) const
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:101
void setIsOnEdge(bool flag, QualWordType &qualWord) const
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:243
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:246
void setTheClu(DetParam const &, ClusterParam &theClusterParam) const
int maxPixelCol() const
virtual LocalPoint localPosition(const MeasurementPoint &) const =0
bool useLAOffsetFromConfig_
Definition: PixelCPEBase.h:231
const TrackerGeometry & geom_
Definition: PixelCPEBase.h:240
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:32
float getLorentzAngle(const uint32_t &) const
float y() const
LocalVector driftDirection
Definition: PixelCPEBase.h:73
const RotationType & rotation() const
bool containsBigPixelInY(int iymin, int iymax) const override
bool containsBigPixelInX(int ixmin, int ixmax) const override
void setHasFilledProb(bool flag, QualWordType &qualWord) const
bool isItEdgePixelInY(int iybin) const override
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
void setQBin(int qbin, QualWordType &qualWord) const
float x() const
bool isTrackerPixel(GeomDetEnumerators::SubDetector m)
float lAWidthFPix_
Definition: PixelCPEBase.h:229
DetParam const & detParam(const GeomDetUnit &det) const