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