CMS 3D CMS Logo

PixelCPEBase.h
Go to the documentation of this file.
1 #ifndef RecoLocalTracker_SiPixelRecHits_interface_PixelCPEBase_h
2 #define RecoLocalTracker_SiPixelRecHits_interface_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 #ifdef EDM_ML_DEBUG
15 #include <atomic>
16 #endif
17 #include <unordered_map>
18 #include <utility>
19 #include <vector>
20 
21 #include <TMath.h>
22 
40 
41 class MagneticField;
43 public:
44  struct DetParam {
45  DetParam() {}
47  // gavril : replace RectangularPixelTopology with PixelTopology
49 
50  //set default value of enum to avoid USBAN errors
53  float theThickness;
54  float thePitchX;
55  float thePitchY;
56 
57  float bz; // local Bz
58  float bx; // local Bx
60  float widthLAFractionX; // Width-LA to Offset-LA in X
61  float widthLAFractionY; // same in Y
62  float lorentzShiftInCmX; // a FULL shift, in cm
63  float lorentzShiftInCmY; // a FULL shift, in cm
64  int detTemplateId; // det if for templates & generic errors
65  int detTemplateId2D; // det if for 2D templates
66  };
67 
68  struct ClusterParam {
71 
72  virtual ~ClusterParam() = default;
73 
74  const SiPixelCluster* theCluster = nullptr;
75 
76  //--- Cluster-level quantities (filled in computeAnglesFrom....)
77  float cotalpha;
78  float cotbeta;
79 
80  // G.Giurgiu (05/14/08) track local coordinates
81  // filled in computeAnglesFrom....
82  float trk_lp_x;
83  float trk_lp_y;
84 
85  // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
86  // with track angles to handle surface deformations (bows/kinks)
87  // filled in computeAnglesFrom.... (btw redundant with the 4 above)
89 
90  //--- Probability (protected by hasFilledProb_)
94  int qBin_; // always filled by qbin
95 
96  bool isOnEdge_; // filled in setTheClu
97  bool hasBadPixels_ = false; // (never used in current code)
98  bool spansTwoROCs_; // filled in setTheClu
99  bool hasFilledProb_ = false;
100  // ggiurgiu@jhu.edu (10/18/2008)
101  bool with_track_angle; // filled in computeAnglesFrom....
102  bool filled_from_2d = false; //
103 
104  // More detailed edge information (for CPE ClusterRepair, and elsewhere...)
105  int edgeTypeX_ = 0; // 0: not on edge, 1: low end on edge, 2: high end
106  int edgeTypeY_ = 0; // 0: not on edge, 1: low end on edge, 2: high end
107  };
108 
109 public:
110  PixelCPEBase(edm::ParameterSet const& conf,
111  const MagneticField* mag,
112  const TrackerGeometry& geom,
113  const TrackerTopology& ttopo,
114  const SiPixelLorentzAngle* lorentzAngle,
115  const SiPixelGenErrorDBObject* genErrorDBObject,
116  const SiPixelTemplateDBObject* templateDBobject,
117  const SiPixelLorentzAngle* lorentzAngleWidth,
118  int flag = 0 // flag=0 for generic, =1 for templates
119  ); // NEW
120 
122 
123  //--------------------------------------------------------------------------
124  // Allow the magnetic field to be set/updated later.
125  //--------------------------------------------------------------------------
126  //inline void setMagField(const MagneticField *mag) const { magfield_ = mag; } // Not used, AH
127 
128  //--------------------------------------------------------------------------
129  // Obtain the angles from the position of the DetUnit.
130  //--------------------------------------------------------------------------
131 
132  inline ReturnType getParameters(const SiPixelCluster& cl, const GeomDetUnit& det) const override {
133 #ifdef EDM_ML_DEBUG
134  nRecHitsTotal_++;
135  LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(all) - " << nRecHitsTotal_;
136 #endif
137 
138  DetParam const& theDetParam = detParam(det);
139  std::unique_ptr<ClusterParam> theClusterParam = createClusterParam(cl);
140  setTheClu(theDetParam, *theClusterParam);
141  computeAnglesFromDetPosition(theDetParam, *theClusterParam);
142 
143  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
144  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
145  LocalError le = localError(theDetParam, *theClusterParam);
146  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
147  auto tuple = std::make_tuple(lp, le, rqw);
148 
149  LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(all) - " << lp.x() << " " << lp.y();
150  return tuple;
151  }
152 
153  //--------------------------------------------------------------------------
154  // In principle we could use the track too to obtain alpha and beta.
155  //--------------------------------------------------------------------------
157  const GeomDetUnit& det,
158  const LocalTrajectoryParameters& ltp) const override {
159 #ifdef EDM_ML_DEBUG
160  nRecHitsTotal_++;
161  LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(on track) - " << nRecHitsTotal_;
162 #endif
163 
164  DetParam const& theDetParam = detParam(det);
165  std::unique_ptr<ClusterParam> theClusterParam = createClusterParam(cl);
166  setTheClu(theDetParam, *theClusterParam);
167  computeAnglesFromTrajectory(theDetParam, *theClusterParam, ltp);
168 
169  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
170  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
171  LocalError le = localError(theDetParam, *theClusterParam);
172  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
173  auto tuple = std::make_tuple(lp, le, rqw);
174 
175  LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(on track) - " << lp.x() << " " << lp.y();
176  return tuple;
177  }
178 
179 private:
180  virtual std::unique_ptr<ClusterParam> createClusterParam(const SiPixelCluster& cl) const = 0;
181 
182  //--------------------------------------------------------------------------
183  // This is where the action happens.
184  //--------------------------------------------------------------------------
185  virtual LocalPoint localPosition(DetParam const& theDetParam, ClusterParam& theClusterParam) const = 0;
186  virtual LocalError localError(DetParam const& theDetParam, ClusterParam& theClusterParam) const = 0;
187 
188  void fillDetParams();
189 
190  //-----------------------------------------------------------------------------
195  //-----------------------------------------------------------------------------
196  SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam& theClusterParam) const;
197 
198 protected:
199  //--- All methods and data members are protected to facilitate (for now)
200  //--- access from derived classes.
201 
203 
204  //---------------------------------------------------------------------------
205  // Data members
206  //---------------------------------------------------------------------------
207 
208  //--- Counters
209 #ifdef EDM_ML_DEBUG
210  mutable std::atomic<int> nRecHitsTotal_; //for debugging only
211  mutable std::atomic<int> nRecHitsUsedEdge_; //for debugging only
212 #endif
213 
214  // Added new members
215  float lAOffset_; // la used to calculate the offset from configuration (for testing)
216  float lAWidthBPix_; // la used to calculate the cluster width from conf.
217  float lAWidthFPix_; // la used to calculate the cluster width from conf.
218  //bool useLAAlignmentOffsets_; // lorentz angle offsets detrmined by alignment
219  bool useLAOffsetFromConfig_; // lorentz angle used to calculate the offset
220  bool useLAWidthFromConfig_; // lorentz angle used to calculate the cluster width
221  bool useLAWidthFromDB_; // lorentz angle used to calculate the cluster width
222 
223  //--- Global quantities
224  int theVerboseLevel; // algorithm's verbosity
225  int theFlag_; // flag to recognice if we are in generic or templates
226 
227  const MagneticField* magfield_; // magnetic field
228  const TrackerGeometry& geom_; // geometry
229  const TrackerTopology& ttopo_; // Tracker Topology
230 
232  const SiPixelLorentzAngle* lorentzAngleWidth_; // for the charge width (generic)
233 
236  bool alpha2Order; // switch on/off E.B effect.
237 
238  bool useLAFromDB_; //Use LA value from the database (used for generic CPE or in template CPE if an error)
239 
242 
243  //errors for template reco for edge hits, based on observed residuals from
244  //studies likely done in 2011...
245  static constexpr float xEdgeXError_ = 23.0f;
246  static constexpr float xEdgeYError_ = 39.0f;
247 
248  static constexpr float yEdgeXError_ = 24.0f;
249  static constexpr float yEdgeYError_ = 96.0f;
250 
251  static constexpr float bothEdgeXError_ = 31.0f;
252  static constexpr float bothEdgeYError_ = 90.0f;
253 
254  static constexpr float clusterSplitMaxError_ = 7777.7f;
255 
256  //---------------------------------------------------------------------------
257  // Geometrical services to subclasses.
258  //---------------------------------------------------------------------------
259 protected:
260  void computeAnglesFromDetPosition(DetParam const& theDetParam, ClusterParam& theClusterParam) const;
261 
262  void computeAnglesFromTrajectory(DetParam const& theDetParam,
263  ClusterParam& theClusterParam,
264  const LocalTrajectoryParameters& ltp) const;
265 
266  void setTheClu(DetParam const&, ClusterParam& theClusterParam) const;
267 
268  LocalVector driftDirection(DetParam& theDetParam, GlobalVector bfield) const;
269  LocalVector driftDirection(DetParam& theDetParam, LocalVector bfield) const;
270  void computeLorentzShifts(DetParam&) const;
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 #endif // RecoLocalTracker_SiPixelRecHits_interface_PixelCPEBase_h
static constexpr float xEdgeYError_
Definition: PixelCPEBase.h:246
Local3DPoint theOrigin
Definition: PixelCPEBase.h:52
void computeAnglesFromTrajectory(DetParam const &theDetParam, ClusterParam &theClusterParam, const LocalTrajectoryParameters &ltp) const
const SiPixelLorentzAngle * lorentzAngleWidth_
Definition: PixelCPEBase.h:232
virtual ~ClusterParam()=default
void computeAnglesFromDetPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:74
ClusterParam(const SiPixelCluster &cl)
Definition: PixelCPEBase.h:70
DetParams m_DetParams
Definition: PixelCPEBase.h:280
static void fillPSetDescription(edm::ParameterSetDescription &desc)
void fillDetParams()
bool doLorentzFromAlignment_
Definition: PixelCPEBase.h:240
static constexpr float bothEdgeYError_
Definition: PixelCPEBase.h:252
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:46
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:51
static constexpr float yEdgeYError_
Definition: PixelCPEBase.h:249
bool useLAWidthFromDB_
Definition: PixelCPEBase.h:221
float lAWidthBPix_
Definition: PixelCPEBase.h:216
void computeLorentzShifts(DetParam &) const
bool useLAWidthFromConfig_
Definition: PixelCPEBase.h:220
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
const PixelTopology * theTopol
Definition: PixelCPEBase.h:48
SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam &theClusterParam) const
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:241
virtual LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const override
Definition: PixelCPEBase.h:132
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:235
const MagneticField * magfield_
Definition: PixelCPEBase.h:227
virtual std::unique_ptr< ClusterParam > createClusterParam(const SiPixelCluster &cl) const =0
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:88
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
static constexpr float xEdgeXError_
Definition: PixelCPEBase.h:245
std::vector< DetParam > DetParams
Definition: PixelCPEBase.h:278
const TrackerTopology & ttopo_
Definition: PixelCPEBase.h:229
const SiPixelLorentzAngle * lorentzAngle_
Definition: PixelCPEBase.h:231
static constexpr float clusterSplitMaxError_
Definition: PixelCPEBase.h:254
const SiPixelGenErrorDBObject * genErrorDBObject_
Definition: PixelCPEBase.h:234
LocalVector driftDirection(DetParam &theDetParam, GlobalVector bfield) const
DetParam const & detParam(const GeomDetUnit &det) const
bool useLAOffsetFromConfig_
Definition: PixelCPEBase.h:219
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:202
Pixel cluster – collection of neighboring pixels above threshold.
const TrackerGeometry & geom_
Definition: PixelCPEBase.h:228
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:27
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det, const LocalTrajectoryParameters &ltp) const override
Definition: PixelCPEBase.h:156
virtual LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
LocalVector driftDirection
Definition: PixelCPEBase.h:59
static constexpr float bothEdgeXError_
Definition: PixelCPEBase.h:251
static constexpr float yEdgeXError_
Definition: PixelCPEBase.h:248
std::tuple< LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType > ReturnType
float lAWidthFPix_
Definition: PixelCPEBase.h:217
void setTheClu(DetParam const &, ClusterParam &theClusterParam) const
#define LogDebug(id)