CMS 3D CMS Logo

PixelCPEBase.h

Go to the documentation of this file.
00001 #ifndef RecoLocalTracker_SiPixelRecHits_PixelCPEBase_H
00002 #define RecoLocalTracker_SiPixelRecHits_PixelCPEBase_H 1
00003 
00004 //-----------------------------------------------------------------------------
00005 // \class        PixelCPEBase
00006 // \description  Base class for pixel CPE's, with all the common code.
00007 // Perform the position and error evaluation of pixel hits using 
00008 // the Det angle to estimate the track impact angle.
00009 // Move geomCorrection to the concrete class. d.k. 06/06.
00010 // change to use Lorentz angle from DB Lotte Wilke, Jan. 31st, 2008
00011 //-----------------------------------------------------------------------------
00012 
00013 #include <utility>
00014 #include <vector>
00015 #include "TMath.h"
00016 
00017 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
00018 #include "RecoLocalTracker/SiPixelRecHits/interface/EtaCorrection.h"
00019 
00020 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00021 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00022 #include "Geometry/TrackerTopology/interface/RectangularPixelTopology.h"
00023 
00024 //--- For the configuration:
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 
00027 //#define TP_OLD
00028 #ifdef TP_OLD
00029 #include "Geometry/CommonDetAlgo/interface/MeasurementPoint.h"
00030 #include "Geometry/CommonDetAlgo/interface/MeasurementError.h"
00031 #include "Geometry/Surface/interface/GloballyPositioned.h"
00032 #else  // new location
00033 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
00034 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementError.h"
00035 #include "DataFormats/GeometrySurface/interface/GloballyPositioned.h"
00036 #endif
00037 
00038 #include "CondFormats/SiPixelObjects/interface/SiPixelLorentzAngle.h"
00039 #include <ext/hash_map>
00040 
00041 class MagneticField;
00042 class SiPixelCPEParmErrors;
00043 class PixelCPEBase : public PixelClusterParameterEstimator {
00044  public:
00045   // PixelCPEBase( const DetUnit& det );
00046   PixelCPEBase(edm::ParameterSet const& conf, const MagneticField * mag = 0, const SiPixelLorentzAngle * lorentzAngle = 0);// const SiPixelCPEParmErrors * parmErrors = 0, const SiPixelLorentzAngle * lorentzAngle = 0);
00047     
00048   //--------------------------------------------------------------------------
00049   // Obtain the angles from the position of the DetUnit.
00050   // LocalValues is typedef for pair<LocalPoint,LocalError> 
00051   //--------------------------------------------------------------------------
00052   inline LocalValues localParameters( const SiPixelCluster & cl, 
00053                                       const GeomDetUnit    & det ) const 
00054   {
00055     nRecHitsTotal_++ ;
00056     setTheDet( det );
00057     computeAnglesFromDetPosition(cl, det);
00058     return std::make_pair( localPosition(cl,det), localError(cl,det) );
00059   }
00060 
00061   //--------------------------------------------------------------------------
00062   // In principle we could use the track too to obtain alpha and beta.
00063   //--------------------------------------------------------------------------
00064   inline LocalValues localParameters( const SiPixelCluster & cl,
00065                                       const GeomDetUnit    & det, 
00066                                       const LocalTrajectoryParameters & ltp) const 
00067   {
00068     nRecHitsTotal_++ ;
00069     setTheDet( det );
00070     computeAnglesFromTrajectory(cl, det, ltp);
00071     return std::make_pair( localPosition(cl,det), localError(cl,det) );
00072   } 
00073 
00074   //--------------------------------------------------------------------------
00075   // The third one, with the user-supplied alpha and beta
00076   //--------------------------------------------------------------------------
00077   inline LocalValues localParameters( const SiPixelCluster & cl,
00078                                       const GeomDetUnit    & det, 
00079                                       float alpha, float beta) const 
00080   {
00081                 nRecHitsTotal_++ ;
00082                 alpha_ = alpha;
00083                 beta_  = beta;
00084                 double HalfPi = 0.5*TMath::Pi();
00085                 cotalpha_ = tan(HalfPi - alpha_);
00086     cotbeta_  = tan(HalfPi - beta_ );
00087                 setTheDet( det );
00088                 return std::make_pair( localPosition(cl,det), localError(cl,det) );
00089   }
00090 
00091 
00092         //--------------------------------------------------------------------------
00093   // Allow the magnetic field to be set/updated later.
00094   //--------------------------------------------------------------------------
00095   inline void setMagField(const MagneticField *mag) const { magfield_ = mag; }
00096 
00097         //--------------------------------------------------------------------------
00098         // Allow the database to be set/updated later.
00099         //--------------------------------------------------------------------------
00100         inline void setDBAccess(const SiPixelCPEParmErrors *parmErrors) const { parmErrors_ = parmErrors; }
00101         
00102         //--------------------------------------------------------------------------
00103   // This is where the action happens.
00104   //--------------------------------------------------------------------------
00105   virtual LocalPoint localPosition(const SiPixelCluster& cl, const GeomDetUnit & det) const;  // = 0, take out dk 8/06
00106   virtual LocalError localError   (const SiPixelCluster& cl, const GeomDetUnit & det) const = 0;
00107   
00108 
00109 
00110   //--------------------------------------------------------------------------
00111   //--- Accessors of other auxiliary quantities
00112   inline float cotAlphaFromCluster() const { return cotAlphaFromCluster_; }
00113   inline float cotBetaFromCluster()  const { return cotBetaFromCluster_; }
00114   inline float probabilityX()        const { return probabilityX_; }
00115   inline float probabilityY()        const { return probabilityY_; }
00116   inline float qBin()                const { return qBin_ ; }
00117 
00118  protected:
00119   //--- All methods and data members are protected to facilitate (for now)
00120   //--- access from derived classes.
00121 
00122   typedef GloballyPositioned<double> Frame;
00123 
00124   //---------------------------------------------------------------------------
00125   //  Data members
00126   //---------------------------------------------------------------------------
00127   //--- Detector-level quantities
00128   mutable const PixelGeomDetUnit * theDet;
00129   mutable const RectangularPixelTopology * theTopol;
00130   mutable GeomDetType::SubDetector thePart;
00131   mutable EtaCorrection theEtaFunc;
00132   mutable float theThickness;
00133   mutable float thePitchX;
00134   mutable float thePitchY;
00135   mutable float theOffsetX;
00136   mutable float theOffsetY;
00137   mutable float theNumOfRow;
00138   mutable float theNumOfCol;
00139   mutable float theDetZ;
00140   mutable float theDetR;
00141   mutable float theLShiftX;
00142   mutable float theLShiftY;
00143   mutable float theSign;
00144 
00145   //--- Cluster-level quantities (may need more)
00146   mutable float alpha_;
00147   mutable float beta_;
00148 
00149   // G.Giurgiu (12/13/06)-----
00150   mutable float cotalpha_;
00151   mutable float cotbeta_;
00152 
00153   // G.Giurgiu (05/14/08) tracl local coordinates
00154   mutable float trk_lp_x;
00155   mutable float trk_lp_y;
00156 
00157 
00158   // [Petar, 5/18/07] 
00159   // Add estimates of cot(alpha) and cot(beta) from the
00160   // cluster length.  This can be used by:
00161   // a) the seed cleaning
00162   // b) any possible crude "quality" flag based on (dis)agreement between
00163   //    W_pred and W (from charge lenght)
00164   // c) an alternative 2nd pass CPE which reads charge per unit length (k_3D) from
00165   //    the DB but then needs angle estimates to switch to 
00166   mutable float cotAlphaFromCluster_;
00167   mutable float cotBetaFromCluster_;
00168 
00169   //--- Probability
00170   mutable float probabilityX_ ; 
00171   mutable float probabilityY_ ; 
00172   mutable float qBin_ ;
00173 
00174   //---------------------------
00175 
00176   // [Petar, 2/23/07]
00177   // Since the sign of the Lorentz shift appears to
00178   // be computed *incorrectly* (i.e. there's a bug) we add new variables
00179   // so that we can study the effect of the bug.
00180   mutable LocalVector driftDirection_;  // drift direction cached // &&&
00181   mutable double lorentzShiftX_;   // a FULL shift, not 1/2 like theLShiftX!
00182   mutable double lorentzShiftY_;   // a FULL shift, not 1/2 like theLShiftY!
00183   mutable double lorentzShiftInCmX_;   // a FULL shift, in cm
00184   mutable double lorentzShiftInCmY_;   // a FULL shift, in cm
00185 
00186   //--- Counters
00187   mutable int    nRecHitsTotal_ ;
00188   mutable int    nRecHitsUsedEdge_ ;
00189 
00190 
00191   //--- Global quantities
00192 //   mutable float theTanLorentzAnglePerTesla;   // tan(Lorentz angle)/Tesla
00193   int     theVerboseLevel;                    // algorithm's verbosity
00194 
00195   mutable const   MagneticField * magfield_;          // magnetic field
00196 
00197         mutable const   SiPixelCPEParmErrors * parmErrors_;
00198         
00199         mutable const SiPixelLorentzAngle * lorentzAngle_;
00200   
00201   bool  alpha2Order;                          // switch on/off E.B effect.
00202 
00203 
00204   //---------------------------------------------------------------------------
00205   //  Methods.
00206   //---------------------------------------------------------------------------
00207   void       setTheDet( const GeomDetUnit & det ) const ;
00208   //
00209   MeasurementPoint measurementPosition( const SiPixelCluster&, 
00210                                         const GeomDetUnit & det) const ;
00211   MeasurementError measurementError   ( const SiPixelCluster&, 
00212                                         const GeomDetUnit & det) const ;
00213 
00214   //---------------------------------------------------------------------------
00215   //  Geometrical services to subclasses.
00216   //---------------------------------------------------------------------------
00217 
00218   //--- Estimation of alpha_ and beta_
00219   //float estimatedAlphaForBarrel(float centerx) const;
00220   void computeAnglesFromDetPosition(const SiPixelCluster & cl, 
00221                                     const GeomDetUnit    & det ) const;
00222   void computeAnglesFromTrajectory (const SiPixelCluster & cl,
00223                                     const GeomDetUnit    & det, 
00224                                     const LocalTrajectoryParameters & ltp) const;
00225   LocalVector driftDirection       ( GlobalVector bfield ) const ; //wrong sign
00226   LocalVector driftDirectionCorrect( GlobalVector bfield ) const ;
00227   void computeLorentzShifts() const ;
00228 
00229   bool isFlipped() const;              // is the det flipped or not?
00230 
00231   //---------------------------------------------------------------------------
00232   //  Cluster-level services.
00233   //---------------------------------------------------------------------------
00234    
00235   //--- Charge on the first, last  and  inner pixels on x and y 
00236   void xCharge(const std::vector<SiPixelCluster::Pixel>&, 
00237                const int&, const int&, float& q1, float& q2) const; 
00238   void yCharge(const std::vector<SiPixelCluster::Pixel>&, 
00239                const int&, const int&, float& q1, float& q2) const; 
00240 
00241   // Temporary fix for older classes
00242   //std::vector<float> 
00243   //xCharge(const std::vector<SiPixelCluster::Pixel>& pixelsVec, 
00244   //    const float& xmin, const float& xmax) const {
00245   // return xCharge(pixelsVec, int(xmin), int(xmax)); 
00246   //}
00247   //std::vector<float> 
00248   // yCharge(const std::vector<SiPixelCluster::Pixel>& pixelsVec, 
00249   //    const float& xmin, const float& xmax) const {
00250   // return yCharge(pixelsVec, int(xmin), int(xmax)); 
00251   //}
00252 
00253 
00254 
00255   //---------------------------------------------------------------------------
00256   //  Various position corrections.
00257   //---------------------------------------------------------------------------
00258   //float geomCorrection()const;
00259 
00260   //--- The Lorentz shift correction
00261   virtual float lorentzShiftX() const;
00262   virtual float lorentzShiftY() const;
00263  
00264   //--- Position in X and Y
00265   virtual float xpos( const SiPixelCluster& ) const = 0;
00266   virtual float ypos( const SiPixelCluster& ) const = 0;
00267   
00268   
00269   
00270  public:
00271   struct Param 
00272   {
00273     Param() : topology(0), drift(0.0, 0.0, 0.0) {}
00274     RectangularPixelTopology const * topology;
00275     LocalVector drift;
00276   };
00277   
00278  private:
00279    typedef  __gnu_cxx::hash_map< unsigned int, Param> Params;
00280   
00281   Params m_Params;
00282   
00283 
00284 
00285 };
00286 
00287 #endif
00288 
00289 

Generated on Tue Jun 9 17:43:59 2009 for CMSSW by  doxygen 1.5.4