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