File indexing completed on 2023-03-17 11:19:38
0001 #ifndef RecoLocalTracker_SiPixelRecHits_interface_PixelCPEBase_h
0002 #define RecoLocalTracker_SiPixelRecHits_interface_PixelCPEBase_h 1
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #ifdef EDM_ML_DEBUG
0015 #include <atomic>
0016 #endif
0017 #include <unordered_map>
0018 #include <utility>
0019 #include <vector>
0020
0021 #include <TMath.h>
0022
0023 #include "CondFormats/SiPixelObjects/interface/SiPixelGenErrorDBObject.h"
0024 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
0025 #include "CondFormats/SiPixelObjects/interface/SiPixelTemplateDBObject.h"
0026 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementError.h"
0027 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
0028 #include "DataFormats/GeometrySurface/interface/GloballyPositioned.h"
0029 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0030 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitQuality.h"
0031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0032 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0035 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0036 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0037 #include "Geometry/CommonTopologies/interface/Topology.h"
0038 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0039 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
0040
0041 class RectangularPixelTopology;
0042 class MagneticField;
0043 class PixelCPEBase : public PixelClusterParameterEstimator {
0044 public:
0045 struct DetParam {
0046 DetParam() {}
0047 const PixelGeomDetUnit* theDet;
0048
0049 const PixelTopology* theTopol;
0050 const RectangularPixelTopology* theRecTopol;
0051
0052
0053 GeomDetType::SubDetector thePart = GeomDetEnumerators::invalidDet;
0054 Local3DPoint theOrigin;
0055 float theThickness;
0056 float thePitchX;
0057 float thePitchY;
0058
0059 float bz;
0060 float bx;
0061 LocalVector driftDirection;
0062 float widthLAFractionX;
0063 float widthLAFractionY;
0064 float lorentzShiftInCmX;
0065 float lorentzShiftInCmY;
0066 int detTemplateId;
0067 int detTemplateId2D;
0068 };
0069
0070 struct ClusterParam {
0071 ClusterParam() {}
0072 ClusterParam(const SiPixelCluster& cl) : theCluster(&cl) {}
0073
0074 virtual ~ClusterParam() = default;
0075
0076 const SiPixelCluster* theCluster = nullptr;
0077
0078
0079 float cotalpha;
0080 float cotbeta;
0081
0082
0083
0084 float trk_lp_x;
0085 float trk_lp_y;
0086
0087
0088
0089
0090 Topology::LocalTrackPred loc_trk_pred;
0091
0092
0093 float probabilityX_;
0094 float probabilityY_;
0095 float probabilityQ_;
0096 int qBin_;
0097
0098 bool isOnEdge_;
0099 bool hasBadPixels_ = false;
0100 bool spansTwoROCs_;
0101 bool hasFilledProb_ = false;
0102
0103 bool with_track_angle;
0104 bool filled_from_2d = false;
0105
0106
0107 int edgeTypeX_ = 0;
0108 int edgeTypeY_ = 0;
0109 };
0110
0111 public:
0112 PixelCPEBase(edm::ParameterSet const& conf,
0113 const MagneticField* mag,
0114 const TrackerGeometry& geom,
0115 const TrackerTopology& ttopo,
0116 const SiPixelLorentzAngle* lorentzAngle,
0117 const SiPixelGenErrorDBObject* genErrorDBObject,
0118 const SiPixelTemplateDBObject* templateDBobject,
0119 const SiPixelLorentzAngle* lorentzAngleWidth,
0120 int flag = 0
0121 );
0122
0123 static void fillPSetDescription(edm::ParameterSetDescription& desc);
0124
0125
0126
0127
0128
0129
0130
0131
0132
0133
0134 inline ReturnType getParameters(const SiPixelCluster& cl, const GeomDetUnit& det) const override {
0135 #ifdef EDM_ML_DEBUG
0136 nRecHitsTotal_++;
0137 LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(all) - " << nRecHitsTotal_;
0138 #endif
0139
0140 DetParam const& theDetParam = detParam(det);
0141 std::unique_ptr<ClusterParam> theClusterParam = createClusterParam(cl);
0142 setTheClu(theDetParam, *theClusterParam);
0143 computeAnglesFromDetPosition(theDetParam, *theClusterParam);
0144
0145
0146 LocalPoint lp = localPosition(theDetParam, *theClusterParam);
0147 LocalError le = localError(theDetParam, *theClusterParam);
0148 SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
0149 auto tuple = std::make_tuple(lp, le, rqw);
0150
0151 LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(all) - " << lp.x() << " " << lp.y();
0152 return tuple;
0153 }
0154
0155
0156
0157
0158 inline ReturnType getParameters(const SiPixelCluster& cl,
0159 const GeomDetUnit& det,
0160 const LocalTrajectoryParameters& ltp) const override {
0161 #ifdef EDM_ML_DEBUG
0162 nRecHitsTotal_++;
0163 LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(on track) - " << nRecHitsTotal_;
0164 #endif
0165
0166 DetParam const& theDetParam = detParam(det);
0167 std::unique_ptr<ClusterParam> theClusterParam = createClusterParam(cl);
0168 setTheClu(theDetParam, *theClusterParam);
0169 computeAnglesFromTrajectory(theDetParam, *theClusterParam, ltp);
0170
0171
0172 LocalPoint lp = localPosition(theDetParam, *theClusterParam);
0173 LocalError le = localError(theDetParam, *theClusterParam);
0174 SiPixelRecHitQuality::QualWordType rqw = rawQualityWord(*theClusterParam);
0175 auto tuple = std::make_tuple(lp, le, rqw);
0176
0177 LogDebug("PixelCPEBase") << " in PixelCPEBase:localParameters(on track) - " << lp.x() << " " << lp.y();
0178 return tuple;
0179 }
0180
0181 private:
0182 virtual std::unique_ptr<ClusterParam> createClusterParam(const SiPixelCluster& cl) const = 0;
0183
0184
0185
0186
0187 virtual LocalPoint localPosition(DetParam const& theDetParam, ClusterParam& theClusterParam) const = 0;
0188 virtual LocalError localError(DetParam const& theDetParam, ClusterParam& theClusterParam) const = 0;
0189
0190 void fillDetParams();
0191
0192
0193
0194
0195
0196
0197
0198 SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam& theClusterParam) const;
0199
0200 protected:
0201
0202
0203
0204 typedef GloballyPositioned<double> Frame;
0205
0206
0207
0208
0209
0210
0211 #ifdef EDM_ML_DEBUG
0212 mutable std::atomic<int> nRecHitsTotal_;
0213 mutable std::atomic<int> nRecHitsUsedEdge_;
0214 #endif
0215
0216
0217 float lAOffset_;
0218 float lAWidthBPix_;
0219 float lAWidthFPix_;
0220
0221 bool useLAOffsetFromConfig_;
0222 bool useLAWidthFromConfig_;
0223 bool useLAWidthFromDB_;
0224
0225
0226 int theVerboseLevel;
0227 int theFlag_;
0228
0229 const MagneticField* magfield_;
0230 const TrackerGeometry& geom_;
0231 const TrackerTopology& ttopo_;
0232
0233 const SiPixelLorentzAngle* lorentzAngle_;
0234 const SiPixelLorentzAngle* lorentzAngleWidth_;
0235
0236 const SiPixelGenErrorDBObject* genErrorDBObject_;
0237 const SiPixelTemplateDBObject* templateDBobject_;
0238 bool alpha2Order;
0239
0240 bool useLAFromDB_;
0241
0242 bool doLorentzFromAlignment_;
0243 bool LoadTemplatesFromDB_;
0244
0245
0246
0247 static constexpr float xEdgeXError_ = 23.0f;
0248 static constexpr float xEdgeYError_ = 39.0f;
0249
0250 static constexpr float yEdgeXError_ = 24.0f;
0251 static constexpr float yEdgeYError_ = 96.0f;
0252
0253 static constexpr float bothEdgeXError_ = 31.0f;
0254 static constexpr float bothEdgeYError_ = 90.0f;
0255
0256 static constexpr float clusterSplitMaxError_ = 7777.7f;
0257
0258
0259
0260
0261 protected:
0262 void computeAnglesFromDetPosition(DetParam const& theDetParam, ClusterParam& theClusterParam) const;
0263
0264 void computeAnglesFromTrajectory(DetParam const& theDetParam,
0265 ClusterParam& theClusterParam,
0266 const LocalTrajectoryParameters& ltp) const;
0267
0268 void setTheClu(DetParam const&, ClusterParam& theClusterParam) const;
0269
0270 LocalVector driftDirection(DetParam& theDetParam, GlobalVector bfield) const;
0271 LocalVector driftDirection(DetParam& theDetParam, LocalVector bfield) const;
0272 void computeLorentzShifts(DetParam&) const;
0273
0274
0275
0276
0277
0278 DetParam const& detParam(const GeomDetUnit& det) const;
0279
0280 using DetParams = std::vector<DetParam>;
0281
0282 DetParams m_DetParams = DetParams(1440);
0283 };
0284
0285 #endif