CMS 3D CMS Logo

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 
27 
28 //--- For the configuration:
30 
31 
35 
37 
38 // new errors
40 // old errors
41 //#include "CondFormats/SiPixelObjects/interface/SiPixelCPEGenericErrorParm.h"
42 
44 
45 #include <unordered_map>
46 
47 #include <iostream>
48 #ifdef EDM_ML_DEBUG
49 #include <atomic>
50 #endif
51 
53 class MagneticField;
55 {
56 public:
57  struct DetParam
58  {
59  DetParam() {}
61  // gavril : replace RectangularPixelTopology with PixelTopology
64 
67  float theThickness;
68  float thePitchX;
69  float thePitchY;
70 
71  float bz; // local Bz
72  float bx; // local Bx
74  float widthLAFractionX; // Width-LA to Offset-LA in X
75  float widthLAFractionY; // same in Y
76  float lorentzShiftInCmX; // a FULL shift, in cm
77  float lorentzShiftInCmY; // a FULL shift, in cm
78  int detTemplateId; // det if for templates & generic errors
79  };
80 
81  struct ClusterParam
82  {
83  ClusterParam(const SiPixelCluster & cl) : theCluster(&cl) {}
84 
85  virtual ~ClusterParam() = default;
86 
88 
89  //--- Cluster-level quantities (filled in computeAnglesFrom....)
90  float cotalpha;
91  float cotbeta;
92 
93  // G.Giurgiu (05/14/08) track local coordinates
94  // filled in computeAnglesFrom....
95  float trk_lp_x;
96  float trk_lp_y;
97 
98  // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
99  // with track angles to handle surface deformations (bows/kinks)
100  // filled in computeAnglesFrom.... (btw redundant with the 4 above)
102 
103  //--- Probability (protected by hasFilledProb_)
107  int qBin_ ; // always filled by qbin
108 
109  bool isOnEdge_ ; // filled in setTheClu
110  bool hasBadPixels_ = false; // (never used in current code)
111  bool spansTwoROCs_ ; // filled in setTheClu
112  bool hasFilledProb_ =false;
113  // ggiurgiu@jhu.edu (10/18/2008)
114  bool with_track_angle; // filled in computeAnglesFrom....
115 
116  };
117 
118 public:
119  PixelCPEBase(edm::ParameterSet const& conf, const MagneticField * mag, const TrackerGeometry& geom, const TrackerTopology& ttopo,
120  const SiPixelLorentzAngle * lorentzAngle,
121  const SiPixelGenErrorDBObject * genErrorDBObject,
122  const SiPixelTemplateDBObject * templateDBobject,
123  const SiPixelLorentzAngle * lorentzAngleWidth,
124  int flag=0 // flag=0 for generic, =1 for templates
125  ); // NEW
126 
127  //--------------------------------------------------------------------------
128  // Allow the magnetic field to be set/updated later.
129  //--------------------------------------------------------------------------
130  //inline void setMagField(const MagneticField *mag) const { magfield_ = mag; } // Not used, AH
131 
132 
133  //--------------------------------------------------------------------------
134  // Obtain the angles from the position of the DetUnit.
135  //--------------------------------------------------------------------------
136 
138  const GeomDetUnit & det ) const
139  {
140 #ifdef EDM_ML_DEBUG
141  nRecHitsTotal_++ ;
142  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<nRecHitsTotal_<<std::endl; //dk
143 #endif
144 
145  DetParam const & theDetParam = detParam(det);
146  ClusterParam * theClusterParam = createClusterParam(cl);
147  setTheClu( theDetParam, *theClusterParam );
148  computeAnglesFromDetPosition(theDetParam, *theClusterParam);
149 
150  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
151  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
152  LocalError le = localError(theDetParam, *theClusterParam);
153  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
154  auto tuple = std::make_tuple(lp, le , rqw);
155  delete theClusterParam;
156 
157  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
158  return tuple;
159  }
160 
161  //--------------------------------------------------------------------------
162  // In principle we could use the track too to obtain alpha and beta.
163  //--------------------------------------------------------------------------
165  const GeomDetUnit & det,
166  const LocalTrajectoryParameters & ltp ) const
167  {
168 #ifdef EDM_ML_DEBUG
169  nRecHitsTotal_++ ;
170  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<nRecHitsTotal_<<std::endl; //dk
171 #endif
172 
173  DetParam const & theDetParam = detParam(det);
174  ClusterParam * theClusterParam = createClusterParam(cl);
175  setTheClu( theDetParam, *theClusterParam );
176  computeAnglesFromTrajectory(theDetParam, *theClusterParam, ltp);
177 
178  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
179  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
180  LocalError le = localError(theDetParam, *theClusterParam);
181  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
182  auto tuple = std::make_tuple(lp, le , rqw);
183  delete theClusterParam;
184 
185  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
186  return tuple;
187  }
188 
189 
190 
191 private:
192  virtual ClusterParam * createClusterParam(const SiPixelCluster & cl) const = 0;
193 
194  //--------------------------------------------------------------------------
195  // This is where the action happens.
196  //--------------------------------------------------------------------------
197  virtual LocalPoint localPosition(DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
198  virtual LocalError localError (DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
199 
200  void fillDetParams();
201 
202  //-----------------------------------------------------------------------------
207  //-----------------------------------------------------------------------------
209 
210 protected:
211  //--- All methods and data members are protected to facilitate (for now)
212  //--- access from derived classes.
213 
215 
216  //---------------------------------------------------------------------------
217  // Data members
218  //---------------------------------------------------------------------------
219 
220  //--- Counters
221 #ifdef EDM_ML_DEBUG
222  mutable std::atomic<int> nRecHitsTotal_ ; //for debugging only
223  mutable std::atomic<int> nRecHitsUsedEdge_ ; //for debugging only
224 #endif
225 
226  // Added new members
227  float lAOffset_; // la used to calculate the offset from configuration (for testing)
228  float lAWidthBPix_; // la used to calculate the cluster width from conf.
229  float lAWidthFPix_; // la used to calculate the cluster width from conf.
230  //bool useLAAlignmentOffsets_; // lorentz angle offsets detrmined by alignment
231  bool useLAOffsetFromConfig_; // lorentz angle used to calculate the offset
232  bool useLAWidthFromConfig_; // lorentz angle used to calculate the cluster width
233  bool useLAWidthFromDB_; // lorentz angle used to calculate the cluster width
234 
235  //--- Global quantities
236  int theVerboseLevel; // algorithm's verbosity
237  int theFlag_; // flag to recognice if we are in generic or templates
238 
239  const MagneticField * magfield_; // magnetic field
240  const TrackerGeometry & geom_; // geometry
241  const TrackerTopology & ttopo_; // Tracker Topology
242 
244  const SiPixelLorentzAngle * lorentzAngleWidth_; // for the charge width (generic)
245 
247  //const SiPixelCPEGenericErrorParm * genErrorParm_; // OLD
248 
250  bool alpha2Order; // switch on/off E.B effect.
251 
254 
255  //---------------------------------------------------------------------------
256  // Geometrical services to subclasses.
257  //---------------------------------------------------------------------------
258 private:
259  void computeAnglesFromDetPosition( DetParam const & theDetParam, ClusterParam & theClusterParam ) const;
260 
261  void computeAnglesFromTrajectory ( DetParam const & theDetParam, ClusterParam & theClusterParam,
262  const LocalTrajectoryParameters & ltp) const;
263 
264  void setTheClu( DetParam const &, ClusterParam & theClusterParam ) const ;
265 
266  LocalVector driftDirection (DetParam & theDetParam, GlobalVector bfield ) const ;
267  LocalVector driftDirection (DetParam & theDetParam, LocalVector bfield ) const ;
268  void computeLorentzShifts(DetParam &) const ;
269 
270  bool isFlipped(DetParam const & theDetParam) const; // is the det flipped or not?
271 
272  //---------------------------------------------------------------------------
273  // Cluster-level services.
274  //---------------------------------------------------------------------------
275 
276  DetParam const & detParam(const GeomDetUnit & det) const;
277 
278  using DetParams=std::vector<DetParam>;
279 
281 
282 };
283 
284 #endif
285 
286 
Local3DPoint theOrigin
Definition: PixelCPEBase.h:66
const SiPixelLorentzAngle * lorentzAngleWidth_
Definition: PixelCPEBase.h:244
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:87
ClusterParam(const SiPixelCluster &cl)
Definition: PixelCPEBase.h:83
bool isFlipped(DetParam const &theDetParam) const
DetParams m_DetParams
Definition: PixelCPEBase.h:280
void fillDetParams()
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
bool useLAWidthFromConfig_
Definition: PixelCPEBase.h:232
const PixelTopology * theTopol
Definition: PixelCPEBase.h:62
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:253
void computeAnglesFromTrajectory(DetParam const &theDetParam, ClusterParam &theClusterParam, const LocalTrajectoryParameters &ltp) const
virtual ClusterParam * createClusterParam(const SiPixelCluster &cl) const =0
void computeLorentzShifts(DetParam &) const
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:249
const MagneticField * magfield_
Definition: PixelCPEBase.h:239
void computeAnglesFromDetPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const
SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam &theClusterParam) const
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:101
std::vector< DetParam > DetParams
Definition: PixelCPEBase.h:278
const TrackerTopology & ttopo_
Definition: PixelCPEBase.h:241
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:243
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:246
void setTheClu(DetParam const &, ClusterParam &theClusterParam) const
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const
Definition: PixelCPEBase.h:164
bool useLAOffsetFromConfig_
Definition: PixelCPEBase.h:231
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:214
Pixel cluster – collection of neighboring pixels above threshold.
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
LocalVector driftDirection
Definition: PixelCPEBase.h:73
std::tuple< LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType > ReturnType
virtual LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
virtual LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const
Definition: PixelCPEBase.h:137
float lAWidthFPix_
Definition: PixelCPEBase.h:229
DetParam const & detParam(const GeomDetUnit &det) const