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