CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelCPEBase.h
Go to the documentation of this file.
1 #ifndef RecoLocalTracker_SiPixelRecHits_PixelCPEBase_H
2 #define RecoLocalTracker_SiPixelRecHits_PixelCPEBase_H 1
3 
4 //-----------------------------------------------------------------------------
5 // \class PixelCPEBase
6 // \description Base class for pixel CPE's, with all the common code.
7 // Perform the position and error evaluation of pixel hits using
8 // the Det angle to estimate the track impact angle.
9 // Move geomCorrection to the concrete class. d.k. 06/06.
10 // change to use Lorentz angle from DB Lotte Wilke, Jan. 31st, 2008
11 // Change to use Generic error & Template calibration from DB - D.Fehling 11/08
12 //-----------------------------------------------------------------------------
13 
14 #include <utility>
15 #include <vector>
16 #include "TMath.h"
17 
20 
25 
26 //--- For the configuration:
28 
29 //#define TP_OLD
30 #ifdef TP_OLD
31 #include "Geometry/CommonDetAlgo/interface/MeasurementPoint.h"
32 #include "Geometry/CommonDetAlgo/interface/MeasurementError.h"
33 #include "Geometry/Surface/interface/GloballyPositioned.h"
34 #else // new location
38 #endif
39 
43 
44 
45 
46 #include <ext/hash_map>
47 
48 #include <iostream>
49 
50 
51 
52 class MagneticField;
54 {
55  public:
56  // PixelCPEBase( const DetUnit& det );
58  const SiPixelLorentzAngle * lorentzAngle = 0, const SiPixelCPEGenericErrorParm * genErrorParm = 0,
59  const SiPixelTemplateDBObject * templateDBobject = 0);
60 
61  //--------------------------------------------------------------------------
62  // Obtain the angles from the position of the DetUnit.
63  // LocalValues is typedef for pair<LocalPoint,LocalError>
64  //--------------------------------------------------------------------------
66  const GeomDetUnit & det ) const
67  {
68  nRecHitsTotal_++ ;
69  setTheDet( det, cl );
71 
72  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
73  LocalPoint lp = localPosition( cl, det );
74  LocalError le = localError( cl, det );
75 
76  return std::make_pair( lp, le );
77  }
78 
79  //--------------------------------------------------------------------------
80  // In principle we could use the track too to obtain alpha and beta.
81  //--------------------------------------------------------------------------
83  const GeomDetUnit & det,
84  const LocalTrajectoryParameters & ltp) const
85  {
86  nRecHitsTotal_++ ;
87  setTheDet( det, cl );
88  computeAnglesFromTrajectory(cl, det, ltp);
89 
90  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
91  LocalPoint lp = localPosition( cl, det );
92  LocalError le = localError( cl, det );
93 
94  return std::make_pair( lp, le );
95  }
96 
97  //--------------------------------------------------------------------------
98  // The third one, with the user-supplied alpha and beta
99  //--------------------------------------------------------------------------
101  const GeomDetUnit & det,
102  float alpha, float beta) const
103  {
104  nRecHitsTotal_++ ;
105  alpha_ = alpha;
106  beta_ = beta;
107  double HalfPi = 0.5*TMath::Pi();
108  cotalpha_ = tan(HalfPi - alpha_);
109  cotbeta_ = tan(HalfPi - beta_ );
110  setTheDet( det, cl );
111 
112  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
113  LocalPoint lp = localPosition( cl, det );
114  LocalError le = localError( cl, det );
115 
116  return std::make_pair( lp, le );
117  }
118 
119 
120  //--------------------------------------------------------------------------
121  // Allow the magnetic field to be set/updated later.
122  //--------------------------------------------------------------------------
123  inline void setMagField(const MagneticField *mag) const { magfield_ = mag; }
124 
125  //--------------------------------------------------------------------------
126  // This is where the action happens.
127  //--------------------------------------------------------------------------
128  virtual LocalPoint localPosition(const SiPixelCluster& cl, const GeomDetUnit & det) const; // = 0, take out dk 8/06
129  virtual LocalError localError (const SiPixelCluster& cl, const GeomDetUnit & det) const = 0;
130 
131 
132 
133  //--------------------------------------------------------------------------
134  //--- Accessors of other auxiliary quantities
135  inline float probabilityX() const { return probabilityX_ ; }
136  inline float probabilityY() const { return probabilityY_ ; }
137  inline float probabilityXY() const {
138  if ( probabilityX_ !=0 && probabilityY_ !=0 )
139  {
141  }
142  else
143  return 0;
144  }
145 
146  inline float probabilityQ() const { return probabilityQ_ ; }
147  inline float qBin() const { return qBin_ ; }
148  inline bool isOnEdge() const { return isOnEdge_ ; }
149  inline bool hasBadPixels() const { return hasBadPixels_ ; }
150  inline bool spansTwoRocks() const { return spansTwoROCs_ ; }
151  inline bool hasFilledProb() const { return hasFilledProb_ ; }
152 
153  //--- Flag to control how SiPixelRecHits compute clusterProbability().
154  //--- Note this is set via the configuration file, and it's simply passed
155  //--- to each TSiPixelRecHit.
156  inline unsigned int clusterProbComputationFlag() const
157  {
159  }
160 
161 
162  //-----------------------------------------------------------------------------
167  //-----------------------------------------------------------------------------
169 
170 
171  protected:
172  //--- All methods and data members are protected to facilitate (for now)
173  //--- access from derived classes.
174 
176 
177  //---------------------------------------------------------------------------
178  // Data members
179  //---------------------------------------------------------------------------
180  //--- Detector-level quantities
181  mutable const PixelGeomDetUnit * theDet;
182 
183  // gavril : replace RectangularPixelTopology with PixelTopology
184  //mutable const RectangularPixelTopology * theTopol;
185  mutable const PixelTopology * theTopol;
186 
187 
189  //mutable EtaCorrection theEtaFunc;
190  mutable float theThickness;
191  mutable float thePitchX;
192  mutable float thePitchY;
193  //mutable float theOffsetX;
194  //mutable float theOffsetY;
195  mutable float theNumOfRow;
196  mutable float theNumOfCol;
197  mutable float theDetZ;
198  mutable float theDetR;
199  mutable float theLShiftX;
200  mutable float theLShiftY;
201  mutable float theSign;
202 
203  //--- Cluster-level quantities (may need more)
204  mutable float alpha_;
205  mutable float beta_;
206 
207  // G.Giurgiu (12/13/06)-----
208  mutable float cotalpha_;
209  mutable float cotbeta_;
210 
211  // G.Giurgiu (05/14/08) track local coordinates
212  mutable float trk_lp_x;
213  mutable float trk_lp_y;
214 
215  //--- Counters
216  mutable int nRecHitsTotal_ ;
217  mutable int nRecHitsUsedEdge_ ;
218 
219  // ggiurgiu@jhu.edu (10/18/2008)
220  mutable bool with_track_angle;
221 
222  //--- Probability
223  mutable float probabilityX_ ;
224  mutable float probabilityY_ ;
225  mutable float probabilityQ_ ;
226  mutable float qBin_ ;
227  mutable bool isOnEdge_ ;
228  mutable bool hasBadPixels_ ;
229  mutable bool spansTwoROCs_ ;
230  mutable bool hasFilledProb_ ;
231 
232  //--- A flag that could be used to change the behavior of
233  //--- clusterProbability() in TSiPixelRecHit (the *transient* one).
234  //--- The problem is that the transient hits are made after the CPE runs
235  //--- and they don't get the access to the PSet, so we pass it via the
236  //--- CPE itself...
237  //
239 
240  //---------------------------
241 
242  // [Petar, 2/23/07]
243  // Since the sign of the Lorentz shift appears to
244  // be computed *incorrectly* (i.e. there's a bug) we add new variables
245  // so that we can study the effect of the bug.
246  mutable LocalVector driftDirection_; // drift direction cached // &&&
247  mutable double lorentzShiftX_; // a FULL shift, not 1/2 like theLShiftX!
248  mutable double lorentzShiftY_; // a FULL shift, not 1/2 like theLShiftY!
249  mutable double lorentzShiftInCmX_; // a FULL shift, in cm
250  mutable double lorentzShiftInCmY_; // a FULL shift, in cm
251 
252 
253  //--- Global quantities
254 // mutable float theTanLorentzAnglePerTesla; // tan(Lorentz angle)/Tesla
255  int theVerboseLevel; // algorithm's verbosity
256 
257  mutable const MagneticField * magfield_; // magnetic field
258 
260 
262 
264 
265  bool alpha2Order; // switch on/off E.B effect.
266 
267  // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
268  // with track angles to handle surface deformations (bows/kinks)
269  //mutable Topology::LocalTrackPred* loc_trk_pred;
271 
273 
274  //---------------------------------------------------------------------------
275  // Methods.
276  //---------------------------------------------------------------------------
277  void setTheDet( const GeomDetUnit & det, const SiPixelCluster & cluster ) const ;
278  //
280  const GeomDetUnit & det) const ;
282  const GeomDetUnit & det) const ;
283 
284  //---------------------------------------------------------------------------
285  // Geometrical services to subclasses.
286  //---------------------------------------------------------------------------
287 
288  //--- Estimation of alpha_ and beta_
289  //float estimatedAlphaForBarrel(float centerx) const;
291  const GeomDetUnit & det ) const;
293  const GeomDetUnit & det,
294  const LocalTrajectoryParameters & ltp) const;
295  LocalVector driftDirection ( GlobalVector bfield ) const ; //wrong sign
297  void computeLorentzShifts() const ;
298 
299  bool isFlipped() const; // is the det flipped or not?
300 
301  //---------------------------------------------------------------------------
302  // Cluster-level services.
303  //---------------------------------------------------------------------------
304 
305  //--- Charge on the first, last and inner pixels on x and y
306  void xCharge(const std::vector<SiPixelCluster::Pixel>&,
307  const int&, const int&, float& q1, float& q2) const;
308  void yCharge(const std::vector<SiPixelCluster::Pixel>&,
309  const int&, const int&, float& q1, float& q2) const;
310 
311  // Temporary fix for older classes
312  //std::vector<float>
313  //xCharge(const std::vector<SiPixelCluster::Pixel>& pixelsVec,
314  // const float& xmin, const float& xmax) const {
315  // return xCharge(pixelsVec, int(xmin), int(xmax));
316  //}
317  //std::vector<float>
318  // yCharge(const std::vector<SiPixelCluster::Pixel>& pixelsVec,
319  // const float& xmin, const float& xmax) const {
320  // return yCharge(pixelsVec, int(xmin), int(xmax));
321  //}
322 
323 
324 
325  //---------------------------------------------------------------------------
326  // Various position corrections.
327  //---------------------------------------------------------------------------
328  //float geomCorrection()const;
329 
330  //--- The Lorentz shift correction
331  virtual float lorentzShiftX() const;
332  virtual float lorentzShiftY() const;
333 
334  //--- Position in X and Y
335  virtual float xpos( const SiPixelCluster& ) const = 0;
336  virtual float ypos( const SiPixelCluster& ) const = 0;
337 
338 
339 
340  public:
341  struct Param
342  {
343  Param() : topology(0), drift(0.0, 0.0, 0.0) {}
344 
345  // giurgiu@jhu.edu 12/09/2010: switch to PixelTopology
346  //RectangularPixelTopology const * topology;
348 
350  };
351 
352  private:
353  typedef __gnu_cxx::hash_map< unsigned int, Param> Params;
354 
356 
357 
358 
359 };
360 
361 #endif
362 
363 
double lorentzShiftX_
Definition: PixelCPEBase.h:247
const double beta
const double Pi
Topology::LocalTrackPred loc_trk_pred_
Definition: PixelCPEBase.h:270
PixelCPEBase(edm::ParameterSet const &conf, const MagneticField *mag=0, const SiPixelLorentzAngle *lorentzAngle=0, const SiPixelCPEGenericErrorParm *genErrorParm=0, const SiPixelTemplateDBObject *templateDBobject=0)
Definition: PixelCPEBase.cc:44
float alpha
Definition: AMPTWrapper.h:95
float theNumOfRow
Definition: PixelCPEBase.h:195
bool with_track_angle
Definition: PixelCPEBase.h:220
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void computeAnglesFromDetPosition(const SiPixelCluster &cl, const GeomDetUnit &det) const
void xCharge(const std::vector< SiPixelCluster::Pixel > &, const int &, const int &, float &q1, float &q2) const
float theLShiftY
Definition: PixelCPEBase.h:200
float theThickness
Definition: PixelCPEBase.h:190
unsigned int clusterProbComputationFlag_
Definition: PixelCPEBase.h:238
SiPixelRecHitQuality::QualWordType rawQualityWord() const
void yCharge(const std::vector< SiPixelCluster::Pixel > &, const int &, const int &, float &q1, float &q2) const
MeasurementError measurementError(const SiPixelCluster &, const GeomDetUnit &det) const
float probabilityQ_
Definition: PixelCPEBase.h:225
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: PixelCPEBase.h:82
bool spansTwoROCs_
Definition: PixelCPEBase.h:229
LocalVector driftDirectionCorrect(GlobalVector bfield) const
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:188
std::pair< LocalPoint, LocalError > LocalValues
double q2[4]
Definition: TauolaWrapper.h:88
virtual float ypos(const SiPixelCluster &) const =0
float theNumOfCol
Definition: PixelCPEBase.h:196
double lorentzShiftInCmY_
Definition: PixelCPEBase.h:250
unsigned int clusterProbComputationFlag() const
Definition: PixelCPEBase.h:156
bool hasBadPixels() const
Definition: PixelCPEBase.h:149
float probabilityY_
Definition: PixelCPEBase.h:224
bool hasBadPixels_
Definition: PixelCPEBase.h:228
LocalVector driftDirection_
Definition: PixelCPEBase.h:246
void computeLorentzShifts() const
bool isFlipped() const
void computeAnglesFromTrajectory(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
bool hasFilledProb() const
Definition: PixelCPEBase.h:151
bool hasFilledProb_
Definition: PixelCPEBase.h:230
LocalVector driftDirection(GlobalVector bfield) const
float probabilityQ() const
Definition: PixelCPEBase.h:146
double lorentzShiftInCmX_
Definition: PixelCPEBase.h:249
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:263
__gnu_cxx::hash_map< unsigned int, Param > Params
Definition: PixelCPEBase.h:353
const MagneticField * magfield_
Definition: PixelCPEBase.h:257
float probabilityX() const
Definition: PixelCPEBase.h:135
void setMagField(const MagneticField *mag) const
Definition: PixelCPEBase.h:123
tuple conf
Definition: dbtoconf.py:185
float probabilityY() const
Definition: PixelCPEBase.h:136
virtual LocalError localError(const SiPixelCluster &cl, const GeomDetUnit &det) const =0
float cl
Definition: Combine.cc:71
LocalTrajectoryParameters loc_traj_param_
Definition: PixelCPEBase.h:272
double q1[4]
Definition: TauolaWrapper.h:87
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:259
const PixelTopology * theTopol
Definition: PixelCPEBase.h:185
MeasurementPoint measurementPosition(const SiPixelCluster &, const GeomDetUnit &det) const
bool isOnEdge() const
Definition: PixelCPEBase.h:148
void setTheDet(const GeomDetUnit &det, const SiPixelCluster &cluster) const
Definition: PixelCPEBase.cc:95
virtual float lorentzShiftX() const
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:175
Pixel cluster – collection of neighboring pixels above threshold.
float probabilityXY() const
Definition: PixelCPEBase.h:137
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const
Definition: PixelCPEBase.h:65
bool spansTwoRocks() const
Definition: PixelCPEBase.h:150
Params m_Params
Definition: PixelCPEBase.h:355
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:181
float theLShiftX
Definition: PixelCPEBase.h:199
PixelTopology const * topology
Definition: PixelCPEBase.h:347
double lorentzShiftY_
Definition: PixelCPEBase.h:248
virtual float lorentzShiftY() const
float probabilityX_
Definition: PixelCPEBase.h:223
virtual LocalPoint localPosition(const SiPixelCluster &cl, const GeomDetUnit &det) const
float qBin() const
Definition: PixelCPEBase.h:147
const SiPixelCPEGenericErrorParm * genErrorParm_
Definition: PixelCPEBase.h:261
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det, float alpha, float beta) const
Definition: PixelCPEBase.h:100
virtual float xpos(const SiPixelCluster &) const =0
int nRecHitsUsedEdge_
Definition: PixelCPEBase.h:217