Back to home page

Project CMSSW displayed by LXR

 
 

    


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 // \class        PixelCPEBase
0006 // \description  Base class for pixel CPE's, with all the common code.
0007 // Perform the position and error evaluation of pixel hits using
0008 // the Det angle to estimate the track impact angle.
0009 // Move geomCorrection to the concrete class. d.k. 06/06.
0010 // change to use Lorentz angle from DB Lotte Wilke, Jan. 31st, 2008
0011 // Change to use Generic error & Template calibration from DB - D.Fehling 11/08
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     // gavril : replace RectangularPixelTopology with PixelTopology
0049     const PixelTopology* theTopol;
0050     const RectangularPixelTopology* theRecTopol;
0051 
0052     //set default value of enum to avoid USBAN errors
0053     GeomDetType::SubDetector thePart = GeomDetEnumerators::invalidDet;
0054     Local3DPoint theOrigin;
0055     float theThickness;
0056     float thePitchX;
0057     float thePitchY;
0058 
0059     float bz;  // local Bz
0060     float bx;  // local Bx
0061     LocalVector driftDirection;
0062     float widthLAFractionX;   // Width-LA to Offset-LA in X
0063     float widthLAFractionY;   // same in Y
0064     float lorentzShiftInCmX;  // a FULL shift, in cm
0065     float lorentzShiftInCmY;  // a FULL shift, in cm
0066     int detTemplateId;        // det if for templates & generic errors
0067     int detTemplateId2D;      // det if for 2D templates
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     //--- Cluster-level quantities (filled in computeAnglesFrom....)
0079     float cotalpha;
0080     float cotbeta;
0081 
0082     // G.Giurgiu (05/14/08) track local coordinates
0083     // filled in computeAnglesFrom....
0084     float trk_lp_x;
0085     float trk_lp_y;
0086 
0087     // ggiurgiu@jhu.edu (12/01/2010) : Needed for calling topology methods
0088     // with track angles to handle surface deformations (bows/kinks)
0089     // filled in computeAnglesFrom.... (btw redundant with the 4 above)
0090     Topology::LocalTrackPred loc_trk_pred;
0091 
0092     //--- Probability  (protected by hasFilledProb_)
0093     float probabilityX_;
0094     float probabilityY_;
0095     float probabilityQ_;
0096     int qBin_;  // always filled by qbin
0097 
0098     bool isOnEdge_;              // filled in setTheClu
0099     bool hasBadPixels_ = false;  // (never used in current code)
0100     bool spansTwoROCs_;          // filled in setTheClu
0101     bool hasFilledProb_ = false;
0102     // ggiurgiu@jhu.edu (10/18/2008)
0103     bool with_track_angle;        // filled in computeAnglesFrom....
0104     bool filled_from_2d = false;  //
0105 
0106     // More detailed edge information (for CPE ClusterRepair, and elsewhere...)
0107     int edgeTypeX_ = 0;  // 0: not on edge, 1: low end on edge, 2: high end
0108     int edgeTypeY_ = 0;  // 0: not on edge, 1: low end on edge, 2: high end
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  // flag=0 for generic, =1 for templates
0121   );                         // NEW
0122 
0123   static void fillPSetDescription(edm::ParameterSetDescription& desc);
0124 
0125   //--------------------------------------------------------------------------
0126   // Allow the magnetic field to be set/updated later.
0127   //--------------------------------------------------------------------------
0128   //inline void setMagField(const MagneticField *mag) const { magfield_ = mag; } // Not used, AH
0129 
0130   //--------------------------------------------------------------------------
0131   // Obtain the angles from the position of the DetUnit.
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     // localPosition( cl, det ) must be called before localError( cl, det ) !!!
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   // In principle we could use the track too to obtain alpha and beta.
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     // localPosition( cl, det ) must be called before localError( cl, det ) !!!
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   // This is where the action happens.
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   //! A convenience method to fill a whole SiPixelRecHitQuality word in one shot.
0194   //! This way, we can keep the details of what is filled within the pixel
0195   //! code and not expose the Transient SiPixelRecHit to it as well.  The name
0196   //! of this function is chosen to match the one in SiPixelRecHit.
0197   //-----------------------------------------------------------------------------
0198   SiPixelRecHitQuality::QualWordType rawQualityWord(ClusterParam& theClusterParam) const;
0199 
0200 protected:
0201   //--- All methods and data members are protected to facilitate (for now)
0202   //--- access from derived classes.
0203 
0204   typedef GloballyPositioned<double> Frame;
0205 
0206   //---------------------------------------------------------------------------
0207   //  Data members
0208   //---------------------------------------------------------------------------
0209 
0210   //--- Counters
0211 #ifdef EDM_ML_DEBUG
0212   mutable std::atomic<int> nRecHitsTotal_;     //for debugging only
0213   mutable std::atomic<int> nRecHitsUsedEdge_;  //for debugging only
0214 #endif
0215 
0216   // Added new members
0217   float lAOffset_;     // la used to calculate the offset from configuration (for testing)
0218   float lAWidthBPix_;  // la used to calculate the cluster width from conf.
0219   float lAWidthFPix_;  // la used to calculate the cluster width from conf.
0220   //bool useLAAlignmentOffsets_; // lorentz angle offsets detrmined by alignment
0221   bool useLAOffsetFromConfig_;  // lorentz angle used to calculate the offset
0222   bool useLAWidthFromConfig_;   // lorentz angle used to calculate the cluster width
0223   bool useLAWidthFromDB_;       // lorentz angle used to calculate the cluster width
0224 
0225   //--- Global quantities
0226   int theVerboseLevel;  // algorithm's verbosity
0227   int theFlag_;         // flag to recognice if we are in generic or templates
0228 
0229   const MagneticField* magfield_;  // magnetic field
0230   const TrackerGeometry& geom_;    // geometry
0231   const TrackerTopology& ttopo_;   // Tracker Topology
0232 
0233   const SiPixelLorentzAngle* lorentzAngle_;
0234   const SiPixelLorentzAngle* lorentzAngleWidth_;  // for the charge width (generic)
0235 
0236   const SiPixelGenErrorDBObject* genErrorDBObject_;
0237   const SiPixelTemplateDBObject* templateDBobject_;
0238   bool alpha2Order;  // switch on/off E.B effect.
0239 
0240   bool useLAFromDB_;  //Use LA value from the database (used for generic CPE or in template CPE if an error)
0241 
0242   bool doLorentzFromAlignment_;
0243   bool LoadTemplatesFromDB_;
0244 
0245   //errors for template reco for edge hits, based on observed residuals from
0246   //studies likely done in 2011...
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   //  Geometrical services to subclasses.
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   //  Cluster-level services.
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  // RecoLocalTracker_SiPixelRecHits_interface_PixelCPEBase_h