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