CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoPixelVertexing/PixelTrackFitting/src/PixelFitterByHelixProjections.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelFitterByHelixProjections.h"
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00007 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00008 
00009 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00010 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00011 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00012 
00013 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00014 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00015 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00016 //#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00017 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00018 
00019 #include "MagneticField/Engine/interface/MagneticField.h"
00020 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00021 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00022 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00023 
00024 #include "CommonTools/Statistics/interface/LinearFit.h"
00025 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00026 
00027 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00028 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00029 #include "FWCore/Framework/interface/ESWatcher.h"
00030 
00031 
00032 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00033 
00034 #include "RZLine.h"
00035 #include "CircleFromThreePoints.h"
00036 #include "PixelTrackBuilder.h"
00037 #include "DataFormats/GeometryVector/interface/Pi.h"
00038 
00039 using namespace std;
00040 
00041 namespace {
00042 
00043   int charge(const std::vector<GlobalPoint> & points) {
00044     // the cross product will tell me...
00045     float dir = (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
00046       - (points[1].y()-points[0].y())*(points[2].x()-points[1].x());
00047     
00048     /*
00049       GlobalVector v21 = points[1]-points[0];
00050       GlobalVector v32 = points[2]-points[1];
00051       float dphi = v32.phi() - v21.phi();
00052       while (dphi >  Geom::fpi()) dphi -=  Geom::ftwoPi();
00053       while (dphi < -Geom::fpi()) dphi +=  Geom::ftwoPi();
00054       return (dphi > 0) ? -1 : 1;
00055     */
00056     return (dir>0) ? -1 : 1;
00057   }
00058 
00059   float cotTheta(const GlobalPoint& inner, const GlobalPoint& outer) {
00060     float dr = outer.perp()-inner.perp();
00061     float dz = outer.z()-inner.z();
00062     return (std::abs(dr) > 1.e-3f) ? dz/dr : 0;
00063   }
00064 
00065   inline float phi(float xC, float yC, int charge) {
00066     return  (charge>0) ? std::atan2(xC,-yC) :  std::atan2(-xC,yC);
00067   }
00068 
00069   float zip(float d0, float phi_p, float curv, 
00070             const GlobalPoint& pinner, const GlobalPoint& pouter) {
00071     //
00072     //phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3) = x(1+x*x/6);
00073     //
00074     
00075     float phi0 = phi_p - Geom::fhalfPi();
00076     GlobalPoint pca(d0*std::cos(phi0), d0*std::sin(phi0),0.);
00077     
00078     float rho2 = curv*curv;
00079     float r1s = (pinner-pca).perp2();
00080     double phi1 = std::sqrt(r1s)*(curv*0.5f)*(1.f+r1s*(rho2/24.f));
00081     float r2s = (pouter-pca).perp2();
00082     double phi2 = std::sqrt(r2s)*(curv*0.5f)*(1.f+r2s*(rho2/24.f));
00083     double z1 = pinner.z();
00084     double z2 = pouter.z();
00085 
00086     if (fabs(curv)>1.e-5) 
00087       return z1 - phi1/(phi1-phi2)*(z1-z2);
00088     else {
00089       double dr = std::max(std::sqrt(r2s)-std::sqrt(r1s),1.e-5f);
00090       return z1-std::sqrt(r1s)*(z2-z1)/dr;
00091     }
00092   }
00093   
00094   double errZip( float apt, float eta) {
00095     double ziperr=0;
00096     double pt = (apt <= 10.) ? apt: 10.;
00097     double p1=0, p2=0,p3=0,p4=0;
00098     float feta = std::abs(eta);
00099     if (feta<=0.8){
00100       p1 = 0.12676e-1;
00101       p2 = -0.22411e-2;
00102       p3 = 0.2987e-3;
00103       p4 = -0.12779e-4;
00104     } else if (feta <=1.6){
00105       p1 = 0.24047e-1;
00106       p2 = -0.66935e-2;
00107       p3 = 0.88111e-3;
00108       p4 = -0.38482e-4;
00109     } else {
00110       p1 = 0.56084e-1;
00111       p2 = -0.13960e-1;
00112       p3 = 0.15744e-2;
00113       p4 = -0.60757e-4;
00114     }
00115     ziperr = p1 + p2*pt + p3*pt*pt +p4*pt*pt*pt;
00116     return ziperr;
00117   }
00118   
00119   double errZip2( float apt, float eta) {
00120     double err = errZip(apt,eta);
00121     return err*err;
00122   }
00123 
00124 
00125   double errTip(float apt, float eta) {
00126     double pt = (apt <= 10.) ? apt : 10.;
00127     double p1=0, p2=0;
00128     float feta = std::abs(eta);
00129     if (feta<=0.8)
00130       {
00131         p1=5.9e-3;
00132         p2=4.7e-3;
00133       }
00134     else if (feta <=1.6){
00135       p1 = 4.9e-3;
00136       p2 = 7.1e-3;
00137     }
00138     else {
00139       p1 = 6.4e-3;
00140       p2 = 1.0e-2;
00141     }
00142     double err=0;
00143     if (pt != 0) err = (p1 + p2/pt);
00144     return err;
00145   }
00146   
00147   double errTip2(float apt, float eta) {
00148     double err = errTip(apt,eta);
00149     return err*err;
00150   }
00151 
00152 
00153 }
00154 
00155 
00156 PixelFitterByHelixProjections::PixelFitterByHelixProjections(
00157    const edm::ParameterSet& cfg) 
00158  : theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0) { }
00159 
00160 reco::Track* PixelFitterByHelixProjections::run(
00161     const edm::EventSetup& es,
00162     const std::vector<const TrackingRecHit * > & hits,
00163     const TrackingRegion & region) const
00164 {
00165   int nhits = hits.size();
00166   if (nhits <2) return 0;
00167 
00168   vector<GlobalPoint> points(nhits);
00169   vector<GlobalError> errors(nhits);
00170   vector<bool> isBarrel(nhits);
00171   
00172   static edm::ESWatcher<TrackerDigiGeometryRecord> watcherTrackerDigiGeometryRecord;
00173   if (!theTracker || watcherTrackerDigiGeometryRecord.check(es)) {
00174     edm::ESHandle<TrackerGeometry> trackerESH;
00175     es.get<TrackerDigiGeometryRecord>().get(trackerESH);
00176     theTracker = trackerESH.product();
00177   }
00178 
00179   static edm::ESWatcher<IdealMagneticFieldRecord>  watcherIdealMagneticFieldRecord;
00180   if (!theField || watcherIdealMagneticFieldRecord.check(es)) {
00181     edm::ESHandle<MagneticField> fieldESH;
00182     es.get<IdealMagneticFieldRecord>().get(fieldESH);
00183     theField = fieldESH.product();
00184   }
00185 
00186   static edm::ESWatcher<TransientRecHitRecord> watcherTransientRecHitRecord;
00187   if (!theTTRecHitBuilder || watcherTransientRecHitRecord.check(es)) {
00188     edm::ESHandle<TransientTrackingRecHitBuilder> ttrhbESH;
00189     std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
00190     es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
00191     theTTRecHitBuilder = ttrhbESH.product();
00192   }
00193 
00194 
00195   for ( int i=0; i!=nhits; ++i) {
00196     TransientTrackingRecHit::RecHitPointer recHit = theTTRecHitBuilder->build(hits[i]);
00197     points[i]  = GlobalPoint( recHit->globalPosition().x()-region.origin().x(), 
00198                               recHit->globalPosition().y()-region.origin().y(),
00199                               recHit->globalPosition().z()-region.origin().z() 
00200                               );
00201     errors[i] = recHit->globalPositionError();
00202     isBarrel[i] = recHit->detUnit()->type().isBarrel();
00203   }
00204 
00205   CircleFromThreePoints circle = (nhits==2) ?
00206         CircleFromThreePoints( GlobalPoint(0.,0.,0.), points[0], points[1]) :
00207         CircleFromThreePoints(points[0],points[1],points[2]); 
00208 
00209   float valPhi, valTip, valPt;
00210 
00211   int iCharge = charge(points);
00212   float curvature = circle.curvature();
00213 
00214   if (curvature > 1.e-4) {
00215     float invPt = PixelRecoUtilities::inversePt( circle.curvature(), es);
00216     valPt = (invPt > 1.e-4f) ? 1.f/invPt : 1.e4f;
00217     CircleFromThreePoints::Vector2D center = circle.center();
00218     valTip = iCharge * (center.mag()-1.f/curvature);
00219     valPhi = phi(center.x(), center.y(), iCharge);
00220   } 
00221   else {
00222     valPt = 1.e4f; 
00223     GlobalVector direction(points[1]-points[0]);
00224     valPhi =  direction.phi(); 
00225     valTip = -points[0].x()*sin(valPhi) + points[0].y()*cos(valPhi); 
00226   }
00227 
00228   float errPt = 0.055f*valPt + 0.017f*valPt*valPt;
00229   float errValTip = errTip(valPt, points.back().eta());
00230   float errPhi = 0.002f;
00231 
00232   float valZip = zip(valTip, valPhi, curvature, points[0],points[1]);
00233   float errValZip = errZip(valPt, points.back().eta());
00234 
00235   float valCotTheta = cotTheta(points[0],points[1]);
00236   float errCotTheta = 0.002f;
00237 
00238   float chi2 = 0;
00239   if (nhits > 2) {
00240     RZLine rzLine(points,errors,isBarrel);
00241     float cottheta, intercept, covss, covii, covsi; 
00242     rzLine.fit(cottheta, intercept, covss, covii, covsi);
00243     chi2 = rzLine.chi2(cottheta, intercept);         //FIXME: check which intercept to use!
00244   }
00245 
00246   PixelTrackBuilder builder;
00247   Measurement1D pt(valPt, errPt);
00248   Measurement1D phi(valPhi, errPhi);
00249   Measurement1D cotTheta(valCotTheta, errCotTheta);
00250   Measurement1D tip(valTip, errValTip);
00251   Measurement1D zip(valZip, errValZip);
00252 
00253   return builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin() );
00254 }
00255 
00256 
00257 
00258