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 
33 
37 
38 
39 
40 #include <unordered_map>
41 
42 #include <iostream>
43 
44 
45 class RectangularPixelTopology;
46 class MagneticField;
48 {
49  public:
50  struct Param
51  {
52  Param() : bz(-9e10f) {}
53  float bz; // local Bz
55  };
56 
57 public:
59  const SiPixelLorentzAngle * lorentzAngle = 0, const SiPixelCPEGenericErrorParm * genErrorParm = 0,
60  const SiPixelTemplateDBObject * templateDBobject = 0);
61 
62  //--------------------------------------------------------------------------
63  // Obtain the angles from the position of the DetUnit.
64  // LocalValues is typedef for pair<LocalPoint,LocalError>
65  //--------------------------------------------------------------------------
67  const GeomDetUnit & det ) const
68  {
69  nRecHitsTotal_++ ;
70  setTheDet( det, cl );
72 
73  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
74  LocalPoint lp = localPosition( cl, det );
75  LocalError le = localError( cl, det );
76 
77  return std::make_pair( lp, le );
78  }
79 
80  //--------------------------------------------------------------------------
81  // In principle we could use the track too to obtain alpha and beta.
82  //--------------------------------------------------------------------------
84  const GeomDetUnit & det,
85  const LocalTrajectoryParameters & ltp) const
86  {
87  nRecHitsTotal_++ ;
88  setTheDet( det, cl );
89  computeAnglesFromTrajectory(cl, det, ltp);
90 
91  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
92  LocalPoint lp = localPosition( cl, det );
93  LocalError le = localError( cl, det );
94 
95  return std::make_pair( lp, le );
96  }
97 
98  //--------------------------------------------------------------------------
99  // The third one, with the user-supplied alpha and beta
100  //--------------------------------------------------------------------------
102  const GeomDetUnit & det,
103  float alpha, float beta) const
104  {
105  nRecHitsTotal_++ ;
106  alpha_ = alpha;
107  beta_ = beta;
108  double HalfPi = 0.5*TMath::Pi();
109  cotalpha_ = tan(HalfPi - alpha_);
110  cotbeta_ = tan(HalfPi - beta_ );
111  setTheDet( det, cl );
112 
113  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
114  LocalPoint lp = localPosition( cl, det );
115  LocalError le = localError( cl, det );
116 
117  return std::make_pair( lp, le );
118  }
119 
120 
122  const GeomDetUnit & det ) const;
123 
124  //--------------------------------------------------------------------------
125  // Allow the magnetic field to be set/updated later.
126  //--------------------------------------------------------------------------
127  inline void setMagField(const MagneticField *mag) const { magfield_ = mag; }
128 
129  //--------------------------------------------------------------------------
130  // This is where the action happens.
131  //--------------------------------------------------------------------------
132  virtual LocalPoint localPosition(const SiPixelCluster& cl, const GeomDetUnit & det) const; // = 0, take out dk 8/06
133  virtual LocalError localError (const SiPixelCluster& cl, const GeomDetUnit & det) const = 0;
134 
135 
136 
137  //--------------------------------------------------------------------------
138  //--- Accessors of other auxiliary quantities
139  inline float probabilityX() const { return probabilityX_ ; }
140  inline float probabilityY() const { return probabilityY_ ; }
141  inline float probabilityXY() const {
142  if ( probabilityX_ !=0 && probabilityY_ !=0 )
143  {
145  }
146  else
147  return 0;
148  }
149 
150  inline float probabilityQ() const { return probabilityQ_ ; }
151  inline float qBin() const { return qBin_ ; }
152  inline bool isOnEdge() const { return isOnEdge_ ; }
153  inline bool hasBadPixels() const { return hasBadPixels_ ; }
154  inline bool spansTwoRocks() const { return spansTwoROCs_ ; }
155  inline bool hasFilledProb() const { return hasFilledProb_ ; }
156 
157  //--- Flag to control how SiPixelRecHits compute clusterProbability().
158  //--- Note this is set via the configuration file, and it's simply passed
159  //--- to each TSiPixelRecHit.
160  inline unsigned int clusterProbComputationFlag() const
161  {
163  }
164 
165 
166  //-----------------------------------------------------------------------------
171  //-----------------------------------------------------------------------------
173 
174 
175  protected:
176  //--- All methods and data members are protected to facilitate (for now)
177  //--- access from derived classes.
178 
180 
181  //---------------------------------------------------------------------------
182  // Data members
183  //---------------------------------------------------------------------------
184  //--- Detector-level quantities
185  mutable const PixelGeomDetUnit * theDet;
186 
187  // gavril : replace RectangularPixelTopology with PixelTopology
188  //mutable const RectangularPixelTopology * theTopol;
189  mutable const PixelTopology * theTopol;
190  mutable const RectangularPixelTopology * theRecTopol;
191 
192  mutable Param const * theParam;
193 
195  //mutable EtaCorrection theEtaFunc;
196  mutable float theThickness;
197  mutable float thePitchX;
198  mutable float thePitchY;
199  //mutable float theOffsetX;
200  //mutable float theOffsetY;
201  mutable float theNumOfRow;
202  mutable float theNumOfCol;
203  mutable float theDetZ;
204  mutable float theDetR;
205  mutable float theLShiftX;
206  mutable float theLShiftY;
207  mutable float theSign;
208 
209  //--- Cluster-level quantities (may need more)
210  mutable float alpha_;
211  mutable float beta_;
212 
213  // G.Giurgiu (12/13/06)-----
214  mutable float cotalpha_;
215  mutable float cotbeta_;
216 
217  // G.Giurgiu (05/14/08) track local coordinates
218  mutable float trk_lp_x;
219  mutable float trk_lp_y;
220 
221  //--- Counters
222  mutable int nRecHitsTotal_ ;
223  mutable int nRecHitsUsedEdge_ ;
224 
225  // ggiurgiu@jhu.edu (10/18/2008)
226  mutable bool with_track_angle;
227 
228  //--- Probability
229  mutable float probabilityX_ ;
230  mutable float probabilityY_ ;
231  mutable float probabilityQ_ ;
232  mutable float qBin_ ;
233  mutable bool isOnEdge_ ;
234  mutable bool hasBadPixels_ ;
235  mutable bool spansTwoROCs_ ;
236  mutable bool hasFilledProb_ ;
237 
238  //--- A flag that could be used to change the behavior of
239  //--- clusterProbability() in TSiPixelRecHit (the *transient* one).
240  //--- The problem is that the transient hits are made after the CPE runs
241  //--- and they don't get the access to the PSet, so we pass it via the
242  //--- CPE itself...
243  //
245 
246  //---------------------------
247 
248  // [Petar, 2/23/07]
249  // Since the sign of the Lorentz shift appears to
250  // be computed *incorrectly* (i.e. there's a bug) we add new variables
251  // so that we can study the effect of the bug.
252  mutable LocalVector driftDirection_; // drift direction cached // &&&
253  mutable double lorentzShiftX_; // a FULL shift, not 1/2 like theLShiftX!
254  mutable double lorentzShiftY_; // a FULL shift, not 1/2 like theLShiftY!
255  mutable double lorentzShiftInCmX_; // a FULL shift, in cm
256  mutable double lorentzShiftInCmY_; // a FULL shift, in cm
257 
258 
259  //--- Global quantities
260 // mutable float theTanLorentzAnglePerTesla; // tan(Lorentz angle)/Tesla
261  int theVerboseLevel; // algorithm's verbosity
262 
263  mutable const MagneticField * magfield_; // magnetic field
264 
266 
268 
270 
271  bool alpha2Order; // switch on/off E.B effect.
272 
273  // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
274  // with track angles to handle surface deformations (bows/kinks)
275  //mutable Topology::LocalTrackPred* loc_trk_pred;
277 
279 
280  //---------------------------------------------------------------------------
281  // Methods.
282  //---------------------------------------------------------------------------
283  void setTheDet( const GeomDetUnit & det, const SiPixelCluster & cluster ) const ;
284 
286  const GeomDetUnit & det) const;
288  const GeomDetUnit & det) const ;
289 
290  //---------------------------------------------------------------------------
291  // Geometrical services to subclasses.
292  //---------------------------------------------------------------------------
293 
295  const GeomDetUnit & det,
296  const LocalTrajectoryParameters & ltp) const;
297  LocalVector driftDirection ( GlobalVector bfield ) const ; //wrong sign
298  LocalVector driftDirection ( LocalVector bfield ) const ; //wrong sign
300  void computeLorentzShifts() const ;
301 
302  bool isFlipped() const; // is the det flipped or not?
303 
304  //---------------------------------------------------------------------------
305  // Cluster-level services.
306  //---------------------------------------------------------------------------
307 
308 
309 
310  //--- The Lorentz shift correction
311  float lorentzShiftX() const;
312  float lorentzShiftY() const;
313 
314  //--- Position in X and Y
315  virtual float xpos( const SiPixelCluster& ) const = 0;
316  virtual float ypos( const SiPixelCluster& ) const = 0;
317 
318 
319  LocalVector const & getDrift() const {return driftDirection_ ;}
320 
321 
322 
323  Param const & param() const;
324 
325  private:
326  typedef std::unordered_map< unsigned int, Param> Params;
327 
328  mutable Params m_Params;
329 
330 
331 
332 };
333 
334 #endif
335 
336 
double lorentzShiftX_
Definition: PixelCPEBase.h:253
const double beta
const double Pi
Topology::LocalTrackPred loc_trk_pred_
Definition: PixelCPEBase.h:276
PixelCPEBase(edm::ParameterSet const &conf, const MagneticField *mag=0, const SiPixelLorentzAngle *lorentzAngle=0, const SiPixelCPEGenericErrorParm *genErrorParm=0, const SiPixelTemplateDBObject *templateDBobject=0)
Definition: PixelCPEBase.cc:45
float alpha
Definition: AMPTWrapper.h:95
float theNumOfRow
Definition: PixelCPEBase.h:201
bool with_track_angle
Definition: PixelCPEBase.h:226
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void computeAnglesFromDetPosition(const SiPixelCluster &cl, const GeomDetUnit &det) const
float theLShiftY
Definition: PixelCPEBase.h:206
float theThickness
Definition: PixelCPEBase.h:196
unsigned int clusterProbComputationFlag_
Definition: PixelCPEBase.h:244
SiPixelRecHitQuality::QualWordType rawQualityWord() const
MeasurementError measurementError(const SiPixelCluster &, const GeomDetUnit &det) const
float probabilityQ_
Definition: PixelCPEBase.h:231
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: PixelCPEBase.h:83
bool spansTwoROCs_
Definition: PixelCPEBase.h:235
LocalVector driftDirectionCorrect(GlobalVector bfield) const
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:194
std::pair< LocalPoint, LocalError > LocalValues
virtual float ypos(const SiPixelCluster &) const =0
float theNumOfCol
Definition: PixelCPEBase.h:202
double lorentzShiftInCmY_
Definition: PixelCPEBase.h:256
unsigned int clusterProbComputationFlag() const
Definition: PixelCPEBase.h:160
bool hasBadPixels() const
Definition: PixelCPEBase.h:153
LocalVector drift
Definition: PixelCPEBase.h:54
float probabilityY_
Definition: PixelCPEBase.h:230
std::unordered_map< unsigned int, Param > Params
Definition: PixelCPEBase.h:326
bool hasBadPixels_
Definition: PixelCPEBase.h:234
LocalVector driftDirection_
Definition: PixelCPEBase.h:252
void computeLorentzShifts() const
bool isFlipped() const
void computeAnglesFromTrajectory(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Param const & param() const
bool hasFilledProb() const
Definition: PixelCPEBase.h:155
bool hasFilledProb_
Definition: PixelCPEBase.h:236
LocalVector driftDirection(GlobalVector bfield) const
float probabilityQ() const
Definition: PixelCPEBase.h:150
double lorentzShiftInCmX_
Definition: PixelCPEBase.h:255
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:269
MeasurementPoint measurementPosition(const SiPixelCluster &cluster, const GeomDetUnit &det) const
const MagneticField * magfield_
Definition: PixelCPEBase.h:263
float probabilityX() const
Definition: PixelCPEBase.h:139
void setMagField(const MagneticField *mag) const
Definition: PixelCPEBase.h:127
tuple conf
Definition: dbtoconf.py:185
float probabilityY() const
Definition: PixelCPEBase.h:140
virtual LocalError localError(const SiPixelCluster &cl, const GeomDetUnit &det) const =0
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:190
LocalTrajectoryParameters loc_traj_param_
Definition: PixelCPEBase.h:278
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:265
const PixelTopology * theTopol
Definition: PixelCPEBase.h:189
bool isOnEdge() const
Definition: PixelCPEBase.h:152
void setTheDet(const GeomDetUnit &det, const SiPixelCluster &cluster) const
Definition: PixelCPEBase.cc:89
float lorentzShiftX() const
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:179
Pixel cluster – collection of neighboring pixels above threshold.
LocalVector const & getDrift() const
Definition: PixelCPEBase.h:319
float probabilityXY() const
Definition: PixelCPEBase.h:141
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const
Definition: PixelCPEBase.h:66
bool spansTwoRocks() const
Definition: PixelCPEBase.h:154
Params m_Params
Definition: PixelCPEBase.h:328
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:185
Param const * theParam
Definition: PixelCPEBase.h:192
float theLShiftX
Definition: PixelCPEBase.h:205
double lorentzShiftY_
Definition: PixelCPEBase.h:254
float lorentzShiftY() const
float probabilityX_
Definition: PixelCPEBase.h:229
virtual LocalPoint localPosition(const SiPixelCluster &cl, const GeomDetUnit &det) const
float qBin() const
Definition: PixelCPEBase.h:151
const SiPixelCPEGenericErrorParm * genErrorParm_
Definition: PixelCPEBase.h:267
LocalValues localParameters(const SiPixelCluster &cl, const GeomDetUnit &det, float alpha, float beta) const
Definition: PixelCPEBase.h:101
virtual float xpos(const SiPixelCluster &) const =0
int nRecHitsUsedEdge_
Definition: PixelCPEBase.h:223