CMS 3D CMS Logo

PixelFitterByConformalMappingAndLine Class Reference

#include <RecoPixelVertexing/PixelTrackFitting/interface/PixelFitterByConformalMappingAndLine.h>

Inheritance diagram for PixelFitterByConformalMappingAndLine:

PixelFitter

List of all members.

Public Member Functions

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

Private Attributes

edm::ParameterSet theConfig


Detailed Description

Definition at line 12 of file PixelFitterByConformalMappingAndLine.h.


Constructor & Destructor Documentation

PixelFitterByConformalMappingAndLine::PixelFitterByConformalMappingAndLine ( const edm::ParameterSet cfg  ) 

Definition at line 40 of file PixelFitterByConformalMappingAndLine.cc.

00041                                 : theConfig(cfg)
00042 { }

PixelFitterByConformalMappingAndLine::PixelFitterByConformalMappingAndLine (  ) 

Definition at line 44 of file PixelFitterByConformalMappingAndLine.cc.

00045 { }

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

Definition at line 16 of file PixelFitterByConformalMappingAndLine.h.

00016 { }


Member Function Documentation

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

Implements PixelFitter.

Definition at line 47 of file PixelFitterByConformalMappingAndLine.cc.

References PixelTrackBuilder::build(), ConformalMappingFit::charge(), RZLine::chi2(), ConformalMappingFit::chi2(), ConformalMappingFit::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), ConformalMappingFit::directionPhi(), e4, Measurement1D::error(), error, HLT_VtxMuL3::errors, edm::ParameterSet::exists(), RZLine::fit(), ConformalMappingFit::fixImpactParmaeter(), edm::EventSetup::get(), edm::ParameterSet::getParameter(), i, ConformalMappingFit::impactParameter(), PixelRecoUtilities::inversePt(), TrackingRegion::origin(), PV3DBase< T, PVType, FrameType >::perp2(), phi, edm::ESHandle< T >::product(), r, GlobalErrorBase< T, ErrorWeightType >::rerr(), sqr(), funct::sqrt(), theConfig, tip, Measurement1D::value(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and z.

00051 {
00052 
00053   int nhits = hits.size();
00054 
00055   vector<GlobalPoint> points;
00056   vector<GlobalError> errors;
00057   vector<bool> isBarrel;
00058   
00059 
00060   edm::ESHandle<TrackerGeometry> tracker;
00061   es.get<TrackerDigiGeometryRecord>().get(tracker);
00062 
00063   edm::ESHandle<MagneticField> field;
00064   es.get<IdealMagneticFieldRecord>().get(field);
00065 
00066   edm::ESHandle<TransientTrackingRecHitBuilder> ttrhBuilder;
00067   string ttrhBuilderName = theConfig.getParameter<std::string>("TTRHBuilder");
00068   es.get<TransientRecHitRecord>().get( ttrhBuilderName, ttrhBuilder);
00069 
00070   for (vector<const TrackingRecHit*>::const_iterator ih=hits.begin();  ih!=hits.end(); ih++) {
00071     TransientTrackingRecHit::RecHitPointer recHit = ttrhBuilder->build(*ih);
00072     points.push_back( recHit->globalPosition() );
00073     errors.push_back( recHit->globalPositionError() );
00074     isBarrel.push_back( recHit->detUnit()->type().isBarrel() );
00075   }
00076 
00077 //    if (useMultScatt) {
00078 //      MultipleScatteringParametrisation ms(hits[i].layer());
00079 //      float cotTheta = (p.z()-zVtx)/p.perp();
00080 //      err += sqr( ms( pt, cotTheta, PixelRecoPointRZ(0.,zVtx) ) );
00081 //    }
00082 
00083   //
00084   // simple fit to get pt, phi0 used for precise calcul.
00085   //
00086   typedef ConformalMappingFit::PointXY PointXY;
00087   vector<PointXY> xy; vector<float> errRPhi2;
00088   for (int i=0; i < nhits; ++i) {
00089     const GlobalPoint & point = points[i];
00090     xy.push_back(PointXY( point.x()-region.origin().x(), point.y()-region.origin().y()));
00091     float phiErr2 = errors[i].phierr(point);
00092     errRPhi2.push_back( point.perp2()*phiErr2);
00093   }
00094   ConformalMappingFit parabola(xy, errRPhi2);
00095   if (theConfig.exists("fixImpactParameter") || nhits < 3)
00096       parabola.fixImpactParmaeter(theConfig.getParameter<double>("fixImpactParameter"));
00097   else if (nhits < 3) {
00098       parabola.fixImpactParmaeter(0.);
00099   }
00100 
00101   Measurement1D curv = parabola.curvature();
00102   float invPt = PixelRecoUtilities::inversePt( curv.value(), es);
00103   float valPt =  (invPt > 1.e-4) ? 1./invPt : 1.e4;
00104   float errPt =PixelRecoUtilities::inversePt(curv.error(), es) * sqr(valPt);
00105   Measurement1D pt (valPt,errPt);
00106   Measurement1D phi = parabola.directionPhi();
00107   Measurement1D tip = parabola.impactParameter();
00108 
00109   //
00110   // precalculate theta to correct errors:
00111   //
00112   vector<float> r(nhits),z(nhits),errZ(nhits);
00113   float simpleCot = ( points.back().z()-points.front().z() )/ (points.back().perp() - points.front().perp() );
00114   for (int i=0; i< nhits; ++i) {
00115     const GlobalPoint & point = points[i]; 
00116     const GlobalError & error = errors[i];
00117     r[i] = sqrt( sqr(point.x()-region.origin().x()) + sqr(point.y()-region.origin().y()) );
00118     r[i] += pixelrecoutilities::LongitudinalBendingCorrection(pt.value(),es)(r[i]);
00119     z[i] = point.z();  
00120     errZ[i] =  (isBarrel[i]) ? sqrt(error.czz()) : sqrt( error.rerr(point) )*simpleCot;
00121   }
00122 
00123   //
00124   // line fit (R-Z plane)
00125   //
00126   RZLine rzLine(r,z,errZ);
00127   float      cottheta, intercept, covss, covii, covsi;
00128   rzLine.fit(cottheta, intercept, covss, covii, covsi);
00129   
00130 //
00131 // parameters for track builder 
00132 //
00133   Measurement1D zip(intercept, sqrt(covii));
00134   Measurement1D cotTheta(cottheta, sqrt(covss));  
00135   float chi2 = parabola.chi2() +  rzLine.chi2(cottheta, intercept);
00136   int charge = parabola.charge();
00137 
00138 
00139   PixelTrackBuilder builder;
00140   return builder.build(pt, phi, cotTheta, tip, zip, chi2, charge, hits,  field.product());
00141 }


Member Data Documentation

edm::ParameterSet PixelFitterByConformalMappingAndLine::theConfig [private]

Definition at line 22 of file PixelFitterByConformalMappingAndLine.h.

Referenced by run().


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