Back to home page

Project CMSSW displayed by LXR

 
 

    


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