CMS 3D CMS Logo

TrackFitter Class Reference

#include <RecoPixelVertexing/PixelLowPtUtilities/interface/TrackFitter.h>

Inheritance diagram for TrackFitter:

PixelFitter

List of all members.

Public Member Functions

virtual reco::Trackrun (const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
 TrackFitter (const edm::ParameterSet &cfg)
virtual ~TrackFitter ()

Private Member Functions

int getCharge (const std::vector< GlobalPoint > &points) const
float getCotThetaAndUpdateZip (const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
void getErrTipAndErrZip (float pt, float eta, float &errZip, float &errTip) const
float getPhi (float xC, float yC, int charge) const
float getZip (float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const

Private Attributes

edm::ParameterSet theConfig
const MagneticFieldtheField
const TrackerGeometrytheTracker
const
TransientTrackingRecHitBuilder
theTTRecHitBuilder


Detailed Description

Definition at line 17 of file TrackFitter.h.


Constructor & Destructor Documentation

TrackFitter::TrackFitter ( const edm::ParameterSet cfg  ) 

Definition at line 40 of file TrackFitter.cc.

00040                                :
00041    theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0)
00042 {
00043 }

virtual TrackFitter::~TrackFitter (  )  [inline, virtual]

Definition at line 21 of file TrackFitter.h.

00021 { }


Member Function Documentation

int TrackFitter::getCharge ( const std::vector< GlobalPoint > &  points  )  const [private]

float TrackFitter::getCotThetaAndUpdateZip ( const GlobalPoint inner,
const GlobalPoint outer,
float  radius,
float  phi,
float  d0,
float &  zip 
) const [private]

Definition at line 153 of file TrackFitter.cc.

References funct::cos(), M_PI_2, muonGeometry::perp(), funct::sin(), and PV3DBase< T, PVType, FrameType >::z().

00155 {
00156    float chi = phi - M_PI_2;
00157    GlobalPoint IP(d0*cos(chi), d0*sin(chi),zip);
00158 
00159    float phi1 = 2*asin(0.5*(inner - IP).perp()/radius);
00160    float phi2 = 2*asin(0.5*(outer - IP).perp()/radius);
00161 
00162    float dr = radius*(phi2 - phi1);
00163    float dz = outer.z()-inner.z();
00164 
00165    // Recalculate ZIP
00166    zip = (inner.z()*phi2 - outer.z()*phi1)/(phi2 - phi1);
00167 
00168    return (fabs(dr) > 1.e-3) ? dz/dr : 0;
00169 }

void TrackFitter::getErrTipAndErrZip ( float  pt,
float  eta,
float &  errZip,
float &  errTip 
) const [private]

Definition at line 205 of file TrackFitter.cc.

References funct::sqrt().

00206 {
00207   float coshEta = cosh(eta);
00208 
00209   { // transverse
00210     float c_ms = 0.0115;
00211     float s_le = 0.0123;
00212     float s_ms2 = c_ms*c_ms / (pt*pt) * coshEta;
00213 
00214     errTip = sqrt(s_le*s_le + s_ms2                  );
00215   }
00216 
00217   { // z
00218     float c_ms = 0.0072;
00219     float s_le = 0.0115;
00220     float s_ms2 = c_ms*c_ms / (pt*pt) * coshEta;
00221 
00222     errZip = sqrt(s_le*s_le + s_ms2 * coshEta*coshEta);
00223   }
00224 }

float TrackFitter::getPhi ( float  xC,
float  yC,
int  charge 
) const [private]

Definition at line 173 of file TrackFitter.cc.

00174 {
00175   float phiC;
00176 
00177   if (charge>0) phiC = atan2(xC,-yC);
00178            else phiC = atan2(-xC,yC);
00179 
00180   return phiC;
00181 }

float TrackFitter::getZip ( float  d0,
float  curv,
const GlobalPoint inner,
const GlobalPoint outer 
) const [private]

Definition at line 185 of file TrackFitter.cc.

References PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::perp2(), r1, r2, and PV3DBase< T, PVType, FrameType >::z().

00187 {
00188   // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
00189   float rho3 = curv*curv*curv;
00190 
00191   float r1 = inner.perp();
00192   double phi1 = r1*curv/2 + inner.perp2()*r1*rho3/48.;
00193 
00194   float r2 = outer.perp();
00195   double phi2 = r2*curv/2 + outer.perp2()*r2*rho3/48.;
00196 
00197   double z1 = inner.z();
00198   double z2 = outer.z();
00199 
00200   return z1 - phi1/(phi1-phi2)*(z1-z2);
00201 }

reco::Track * TrackFitter::run ( const edm::EventSetup es,
const std::vector< const TrackingRecHit * > &  hits,
const TrackingRegion region 
) const [virtual]

Implements PixelFitter.

Definition at line 47 of file TrackFitter.cc.

References TransientTrackingRecHitBuilder::build(), PixelTrackBuilder::build(), CircleFromThreePoints::center(), RZLine::chi2(), PixelRecoUtilities::curvature(), CircleFromThreePoints::curvature(), e4, HLT_VtxMuL3::errors, RZLine::fit(), edm::EventSetup::get(), PixelRecoUtilities::inversePt(), Basic2DVector< T >::mag(), phi, edm::ESHandle< T >::product(), tip, Basic2DVector< T >::x(), and Basic2DVector< T >::y().

00050 {
00051   int nhits = hits.size();
00052   if(nhits <2) return 0;
00053 
00054   vector<GlobalPoint> points;
00055   vector<GlobalError> errors;
00056   vector<bool> isBarrel;
00057   
00058   if (!theField || !theTracker || !theTTRecHitBuilder)
00059   {
00060     edm::ESHandle<TrackerGeometry> trackerESH;
00061     es.get<TrackerDigiGeometryRecord>().get(trackerESH);
00062     theTracker = trackerESH.product();
00063 
00064     edm::ESHandle<MagneticField> fieldESH;
00065     es.get<IdealMagneticFieldRecord>().get(fieldESH);
00066     theField = fieldESH.product();
00067 
00068     edm::ESHandle<TransientTrackingRecHitBuilder> ttrhbESH;
00069     std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
00070     es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
00071     theTTRecHitBuilder = ttrhbESH.product();
00072   }
00073 
00074   for (vector<const TrackingRecHit*>::const_iterator ih = hits.begin();
00075                                                      ih!= hits.end(); ih++)
00076   {
00077     TransientTrackingRecHit::RecHitPointer recHit =
00078       theTTRecHitBuilder->build(*ih);
00079 
00080     points.push_back(recHit->globalPosition());
00081     errors.push_back(recHit->globalPositionError());
00082     isBarrel.push_back(recHit->detUnit()->type().isBarrel() );
00083   }
00084 
00085   CircleFromThreePoints circle = (nhits==2) ?
00086         CircleFromThreePoints(GlobalPoint(0.,0.,0.), points[0], points[1]) :
00087         CircleFromThreePoints(points[0],points[1],points[2]); 
00088 
00089   int charge = getCharge(points);
00090   float curvature = circle.curvature();
00091 
00092   // pt
00093   float invPt = PixelRecoUtilities::inversePt(curvature, es);
00094   float valPt = (invPt > 1.e-4) ? 1./invPt : 1.e4;
00095   float errPt = 0.055*valPt + 0.017*valPt*valPt;
00096 
00097   CircleFromThreePoints::Vector2D center = circle.center();
00098 
00099   // tip
00100   float valTip = charge * (center.mag()-1/curvature);
00101   // zip
00102   float valZip = getZip(valTip, curvature, points[0],points[1]);
00103   // phi
00104   float valPhi = getPhi(center.x(), center.y(), charge);
00105   // cot(theta), update zip
00106   float valCotTheta =
00107      getCotThetaAndUpdateZip(points[0],points[1], 1/curvature,
00108                              valPhi,valTip,valZip);
00109 
00110   // errors
00111   float errTip, errZip;
00112   getErrTipAndErrZip(valPt, points.back().eta(), errTip,errZip);
00113   float errPhi      = 0.002;
00114   float errCotTheta = 0.002;
00115 
00116   float chi2 = 0;
00117   if(nhits > 2)
00118   {
00119     RZLine rzLine(points,errors,isBarrel);
00120     float      cotTheta, intercept, covss, covii, covsi; 
00121     rzLine.fit(cotTheta, intercept, covss, covii, covsi);
00122     chi2 = rzLine.chi2(cotTheta, intercept);
00123     //FIXME: check which intercept to use!
00124   }
00125 
00126   // build pixel track
00127   PixelTrackBuilder builder;
00128 
00129   Measurement1D pt      (valPt,       errPt);
00130   Measurement1D phi     (valPhi,      errPhi);
00131   Measurement1D cotTheta(valCotTheta, errCotTheta);
00132   Measurement1D tip     (valTip,      errTip);
00133   Measurement1D zip     (valZip,      errZip);
00134 
00135   return builder.build(pt, phi, cotTheta, tip, zip, chi2,
00136                        charge, hits, theField);
00137 }


Member Data Documentation

edm::ParameterSet TrackFitter::theConfig [private]

Definition at line 38 of file TrackFitter.h.

const MagneticField* TrackFitter::theField [mutable, private]

Definition at line 41 of file TrackFitter.h.

const TrackerGeometry* TrackFitter::theTracker [mutable, private]

Definition at line 40 of file TrackFitter.h.

const TransientTrackingRecHitBuilder* TrackFitter::theTTRecHitBuilder [mutable, private]

Definition at line 42 of file TrackFitter.h.


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:34:11 2009 for CMSSW by  doxygen 1.5.4