CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
26 
27 //--- For the configuration:
29 
30 
34 
36 
37 // new errors
39 // old errors
40 //#include "CondFormats/SiPixelObjects/interface/SiPixelCPEGenericErrorParm.h"
41 
43 
44 #include <unordered_map>
45 
46 #include <iostream>
47 #ifdef EDM_ML_DEBUG
48 #include <atomic>
49 #endif
50 
52 class MagneticField;
54 {
55  public:
56  struct DetParam
57  {
58  DetParam() {}
60  // gavril : replace RectangularPixelTopology with PixelTopology
63 
66  float theThickness;
67  float thePitchX;
68  float thePitchY;
69  //float theDetR;
70  //float theDetZ;
71  //float theNumOfRow; //Not used, AH
72  //float theNumOfCol; //Not used, AH
73  //float theSign; //Not used, AH
74 
75  float bz; // local Bz
77  float widthLAFractionX; // Width-LA to Offset-LA in X
78  float widthLAFractionY; // same in Y
79  float lorentzShiftInCmX; // a FULL shift, in cm
80  float lorentzShiftInCmY; // a FULL shift, in cm
81  int detTemplateId; // det if for templates & generic errors
82  };
83 
84  struct ClusterParam
85  {
86  ClusterParam(const SiPixelCluster & cl) : theCluster(&cl), loc_trk_pred(0.0,0.0,0.0,0.0),
87  probabilityX_(0.0), probabilityY_(0.0), probabilityQ_(0.0), qBin_(0.0),
90 
91  //--- Cluster-level quantities (may need more)
92  float cotalpha;
93  float cotbeta;
94  //bool zneg; // Not used, AH
95 
96  // G.Giurgiu (05/14/08) track local coordinates
97  float trk_lp_x;
98  float trk_lp_y;
99 
100  // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
101  // with track angles to handle surface deformations (bows/kinks)
103  //LocalTrajectoryParameters loc_traj_param; // Not used, AH
104 
105  // ggiurgiu@jhu.edu (10/18/2008)
107 
108  //--- Probability
109  float probabilityX_ ;
110  float probabilityY_ ;
111  float probabilityQ_ ;
112  float qBin_ ;
113  bool isOnEdge_ ;
117  };
118 
119 public:
121  const SiPixelLorentzAngle * lorentzAngle,
122  const SiPixelGenErrorDBObject * genErrorDBObject,
123  const SiPixelTemplateDBObject * templateDBobject,
124  const SiPixelLorentzAngle * lorentzAngleWidth,
125  int flag=0 // flag=0 for generic, =1 for templates
126  ); // NEW
127 
128  //--------------------------------------------------------------------------
129  // Allow the magnetic field to be set/updated later.
130  //--------------------------------------------------------------------------
131  //inline void setMagField(const MagneticField *mag) const { magfield_ = mag; } // Not used, AH
132 
133 
134  //--------------------------------------------------------------------------
135  // Obtain the angles from the position of the DetUnit.
136  //--------------------------------------------------------------------------
137 
139  const GeomDetUnit & det ) const
140  {
141 #ifdef EDM_ML_DEBUG
142  nRecHitsTotal_++ ;
143  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<nRecHitsTotal_<<std::endl; //dk
144 #endif
145 
146  DetParam const & theDetParam = detParam(det);
147  ClusterParam * theClusterParam = createClusterParam(cl);
148  setTheClu( theDetParam, *theClusterParam );
149  computeAnglesFromDetPosition(theDetParam, *theClusterParam);
150 
151  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
152  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
153  LocalError le = localError(theDetParam, *theClusterParam);
154  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
155  auto tuple = std::make_tuple(lp, le , rqw);
156  delete theClusterParam;
157 
158  //std::cout<<" in PixelCPEBase:localParameters(all) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
159  return tuple;
160  }
161 
162  //--------------------------------------------------------------------------
163  // In principle we could use the track too to obtain alpha and beta.
164  //--------------------------------------------------------------------------
166  const GeomDetUnit & det,
167  const LocalTrajectoryParameters & ltp ) const
168  {
169 #ifdef EDM_ML_DEBUG
170  nRecHitsTotal_++ ;
171  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<nRecHitsTotal_<<std::endl; //dk
172 #endif
173 
174  DetParam const & theDetParam = detParam(det);
175  ClusterParam * theClusterParam = createClusterParam(cl);
176  setTheClu( theDetParam, *theClusterParam );
177  computeAnglesFromTrajectory(theDetParam, *theClusterParam, ltp);
178 
179  // localPosition( cl, det ) must be called before localError( cl, det ) !!!
180  LocalPoint lp = localPosition(theDetParam, *theClusterParam);
181  LocalError le = localError(theDetParam, *theClusterParam);
182  SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
183  auto tuple = std::make_tuple(lp, le , rqw);
184  delete theClusterParam;
185 
186  //std::cout<<" in PixelCPEBase:localParameters(on track) - "<<lp.x()<<" "<<lp.y()<<std::endl; //dk
187  return tuple;
188  }
189 
190 
191 
192 private:
193  virtual ClusterParam * createClusterParam(const SiPixelCluster & cl) const = 0;
194 
195  //--------------------------------------------------------------------------
196  // This is where the action happens.
197  //--------------------------------------------------------------------------
198  virtual LocalPoint localPosition(DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
199  virtual LocalError localError (DetParam const & theDetParam, ClusterParam & theClusterParam) const = 0;
200 
201  void fillDetParams();
202 
203  //-----------------------------------------------------------------------------
208  //-----------------------------------------------------------------------------
209  SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam & theClusterParam) const;
210 
211  protected:
212  //--- All methods and data members are protected to facilitate (for now)
213  //--- access from derived classes.
214 
216 
217  //---------------------------------------------------------------------------
218  // Data members
219  //---------------------------------------------------------------------------
220 
221  //--- Counters
222 #ifdef EDM_ML_DEBUG
223  mutable std::atomic<int> nRecHitsTotal_ ; //for debugging only
224  mutable std::atomic<int> nRecHitsUsedEdge_ ; //for debugging only
225 #endif
226 
227  // Added new members
228  float lAOffset_; // la used to calculate the offset from configuration (for testing)
229  float lAWidthBPix_; // la used to calculate the cluster width from conf.
230  float lAWidthFPix_; // la used to calculate the cluster width from conf.
231  //bool useLAAlignmentOffsets_; // lorentz angle offsets detrmined by alignment
232  bool useLAOffsetFromConfig_; // lorentz angle used to calculate the offset
233  bool useLAWidthFromConfig_; // lorentz angle used to calculate the cluster width
234  bool useLAWidthFromDB_; // lorentz angle used to calculate the cluster width
235 
236  //--- Global quantities
237  int theVerboseLevel; // algorithm's verbosity
238  int theFlag_; // flag to recognice if we are in generic or templates
239 
240  const MagneticField * magfield_; // magnetic field
241  const TrackerGeometry & geom_; // geometry
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:65
const SiPixelLorentzAngle * lorentzAngleWidth_
Definition: PixelCPEBase.h:244
LocalVector driftDirection(DetParam &theDetParam, GlobalVector bfield) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const SiPixelCluster * theCluster
Definition: PixelCPEBase.h:89
ClusterParam(const SiPixelCluster &cl)
Definition: PixelCPEBase.h:86
bool isFlipped(DetParam const &theDetParam) const
DetParams m_DetParams
Definition: PixelCPEBase.h:280
void fillDetParams()
const PixelGeomDetUnit * theDet
Definition: PixelCPEBase.h:59
GeomDetType::SubDetector thePart
Definition: PixelCPEBase.h:64
const RectangularPixelTopology * theRecTopol
Definition: PixelCPEBase.h:62
bool useLAWidthFromDB_
Definition: PixelCPEBase.h:234
float lAWidthBPix_
Definition: PixelCPEBase.h:229
bool useLAWidthFromConfig_
Definition: PixelCPEBase.h:233
virtual ClusterParam * createClusterParam(const SiPixelCluster &cl) const =0
const PixelTopology * theTopol
Definition: PixelCPEBase.h:61
bool LoadTemplatesFromDB_
Definition: PixelCPEBase.h:253
void computeAnglesFromTrajectory(DetParam const &theDetParam, ClusterParam &theClusterParam, const LocalTrajectoryParameters &ltp) const
virtual LocalPoint localPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
void computeLorentzShifts(DetParam &) const
const SiPixelTemplateDBObject * templateDBobject_
Definition: PixelCPEBase.h:249
const MagneticField * magfield_
Definition: PixelCPEBase.h:240
void computeAnglesFromDetPosition(DetParam const &theDetParam, ClusterParam &theClusterParam) const
SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam &theClusterParam) const
tuple conf
Definition: dbtoconf.py:185
Topology::LocalTrackPred loc_trk_pred
Definition: PixelCPEBase.h:102
std::vector< DetParam > DetParams
Definition: PixelCPEBase.h:278
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:165
bool useLAOffsetFromConfig_
Definition: PixelCPEBase.h:232
PixelCPEBase(edm::ParameterSet const &conf, const MagneticField *mag, const TrackerGeometry &geom, const SiPixelLorentzAngle *lorentzAngle, const SiPixelGenErrorDBObject *genErrorDBObject, const SiPixelTemplateDBObject *templateDBobject, const SiPixelLorentzAngle *lorentzAngleWidth, int flag=0)
Definition: PixelCPEBase.cc:44
GloballyPositioned< double > Frame
Definition: PixelCPEBase.h:215
Pixel cluster – collection of neighboring pixels above threshold.
const TrackerGeometry & geom_
Definition: PixelCPEBase.h:241
virtual LocalError localError(DetParam const &theDetParam, ClusterParam &theClusterParam) const =0
LocalVector driftDirection
Definition: PixelCPEBase.h:76
volatile std::atomic< bool > shutdown_flag false
std::tuple< LocalPoint, LocalError, SiPixelRecHitQuality::QualWordType > ReturnType
ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const
Definition: PixelCPEBase.h:138
float lAWidthFPix_
Definition: PixelCPEBase.h:230
DetParam const & detParam(const GeomDetUnit &det) const