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  bool filled_from_2d = false; //
116 
117  // More detailed edge information (for CPE ClusterRepair, and elsewhere...)
118  int edgeTypeX_ = 0; // 0: not on edge, 1: low end on edge, 2: high end
119  int edgeTypeY_ = 0; // 0: not on edge, 1: low end on edge, 2: high end
120  };
121 
122 public:
123  PixelCPEBase(edm::ParameterSet const& conf, const MagneticField * mag, const TrackerGeometry& geom, const TrackerTopology& ttopo,
124  const SiPixelLorentzAngle * lorentzAngle,
125  const SiPixelGenErrorDBObject * genErrorDBObject,
126  const SiPixelTemplateDBObject * templateDBobject,
127  const SiPixelLorentzAngle * lorentzAngleWidth,
128  int flag=0 // flag=0 for generic, =1 for templates
129  ); // NEW
130 
131  //--------------------------------------------------------------------------
132  // Allow the magnetic field to be set/updated later.
133  //--------------------------------------------------------------------------
134  //inline void setMagField(const MagneticField *mag) const { magfield_ = mag; } // Not used, AH
135 
136 
137  //--------------------------------------------------------------------------
138  // Obtain the angles from the position of the DetUnit.
139  //--------------------------------------------------------------------------
140 
142  const GeomDetUnit & det ) const override
143  {
144 #ifdef EDM_ML_DEBUG
145  nRecHitsTotal_++ ;
146  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<nRecHitsTotal_<<std::endl; //dk
147 #endif
148 
149  DetParam const & theDetParam = detParam(det);
150  ClusterParam * theClusterParam = createClusterParam(cl);
151  setTheClu( theDetParam, *theClusterParam );
152  computeAnglesFromDetPosition(theDetParam, *theClusterParam);
153 
154  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
155  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
156  LocalError le = localError(theDetParam, *theClusterParam);
157  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
158  auto tuple = std::make_tuple(lp, le , rqw);
159  delete theClusterParam;
160 
161  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
162  return tuple;
163  }
164 
165  //--------------------------------------------------------------------------
166  // In principle we could use the track too to obtain alpha and beta.
167  //--------------------------------------------------------------------------
169  const GeomDetUnit & det,
170  const LocalTrajectoryParameters & ltp ) const override
171  {
172 #ifdef EDM_ML_DEBUG
173  nRecHitsTotal_++ ;
174  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<nRecHitsTotal_<<std::endl; //dk
175 #endif
176 
177  DetParam const & theDetParam = detParam(det);
178  ClusterParam * theClusterParam = createClusterParam(cl);
179  setTheClu( theDetParam, *theClusterParam );
180  computeAnglesFromTrajectory(theDetParam, *theClusterParam, ltp);
181 
182  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
183  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
184  LocalError le = localError(theDetParam, *theClusterParam);
185  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
186  auto tuple = std::make_tuple(lp, le , rqw);
187  delete theClusterParam;
188 
189  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
190  return tuple;
191  }
192 
193 
194 
195 private:
196  virtual ClusterParam * createClusterParam(const SiPixelCluster & cl) const = 0;
197 
198  //--------------------------------------------------------------------------
199  // This is where the action happens.
200  //--------------------------------------------------------------------------
201  virtual LocalPoint localPosition(DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
202  virtual LocalError localError (DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
203 
204  void fillDetParams();
205 
206  //-----------------------------------------------------------------------------
211  //-----------------------------------------------------------------------------
213 
214 protected:
215  //--- All methods and data members are protected to facilitate (for now)
216  //--- access from derived classes.
217 
219 
220  //---------------------------------------------------------------------------
221  // Data members
222  //---------------------------------------------------------------------------
223 
224  //--- Counters
225 #ifdef EDM_ML_DEBUG
226  mutable std::atomic<int> nRecHitsTotal_ ; //for debugging only
227  mutable std::atomic<int> nRecHitsUsedEdge_ ; //for debugging only
228 #endif
229 
230  // Added new members
231  float lAOffset_; // la used to calculate the offset from configuration (for testing)
232  float lAWidthBPix_; // la used to calculate the cluster width from conf.
233  float lAWidthFPix_; // la used to calculate the cluster width from conf.
234  //bool useLAAlignmentOffsets_; // lorentz angle offsets detrmined by alignment
235  bool useLAOffsetFromConfig_; // lorentz angle used to calculate the offset
236  bool useLAWidthFromConfig_; // lorentz angle used to calculate the cluster width
237  bool useLAWidthFromDB_; // lorentz angle used to calculate the cluster width
238 
239  //--- Global quantities
240  int theVerboseLevel; // algorithm's verbosity
241  int theFlag_; // flag to recognice if we are in generic or templates
242 
243  const MagneticField * magfield_; // magnetic field
244  const TrackerGeometry & geom_; // geometry
245  const TrackerTopology & ttopo_; // Tracker Topology
246 
248  const SiPixelLorentzAngle * lorentzAngleWidth_; // for the charge width (generic)
249 
251  //const SiPixelCPEGenericErrorParm * genErrorParm_; // OLD
252 
254  bool alpha2Order; // switch on/off E.B effect.
255 
258 
259  //---------------------------------------------------------------------------
260  // Geometrical services to subclasses.
261  //---------------------------------------------------------------------------
262 protected:
263  void computeAnglesFromDetPosition( DetParam const & theDetParam, ClusterParam & theClusterParam ) const;
264 
265  void computeAnglesFromTrajectory ( DetParam const & theDetParam, ClusterParam & theClusterParam,
266  const LocalTrajectoryParameters & ltp) const;
267 
268  void setTheClu( DetParam const &, ClusterParam & theClusterParam ) const ;
269 
270  LocalVector driftDirection (DetParam & theDetParam, GlobalVector bfield ) const ;
271  LocalVector driftDirection (DetParam & theDetParam, LocalVector bfield ) const ;
272  void computeLorentzShifts(DetParam &) const ;
273 
274 
275  //---------------------------------------------------------------------------
276  // Cluster-level services.
277  //---------------------------------------------------------------------------
278 
279  DetParam const & detParam(const GeomDetUnit & det) const;
280 
281  using DetParams=std::vector<DetParam>;
282 
284 
285 };
286 
287 #endif
288 
289 
Local3DPoint theOrigin
Definition: PixelCPEBase.h:66
const SiPixelLorentzAngle * lorentzAngleWidth_
Definition: PixelCPEBase.h:248
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const override
Definition: PixelCPEBase.h:141
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
DetParams m_DetParams
Definition: PixelCPEBase.h:283
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:237
float lAWidthBPix_
Definition: PixelCPEBase.h:232
bool useLAWidthFromConfig_
Definition: PixelCPEBase.h:236
const PixelTopology * theTopol
Definition: PixelCPEBase.h:62
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:257
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:253
const MagneticField * magfield_
Definition: PixelCPEBase.h:243
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:281
const TrackerTopology & ttopo_
Definition: PixelCPEBase.h:245
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:247
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:250
void setTheClu(DetParam const &, ClusterParam &theClusterParam) const
bool useLAOffsetFromConfig_
Definition: PixelCPEBase.h:235
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:218
Pixel cluster – collection of neighboring pixels above threshold.
const TrackerGeometry & geom_
Definition: PixelCPEBase.h:244
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
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const override
Definition: PixelCPEBase.h:168
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
float lAWidthFPix_
Definition: PixelCPEBase.h:233
DetParam const & detParam(const GeomDetUnit &det) const