CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoPixelVertexing/PixelTrackFitting/src/KFBasedPixelFitter.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelTrackFitting/interface/KFBasedPixelFitter.h"
00002 
00003 #include "FWCore/Framework/interface/EventSetup.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 
00010 #include "MagneticField/Engine/interface/MagneticField.h"
00011 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00012 
00013 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHitBuilder.h"
00014 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00015 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00016 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00017 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00018 
00019 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00020 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00021 
00022 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00023 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00025 
00026 #include "DataFormats/TrackReco/interface/Track.h"
00027 #include "DataFormats/TrackReco/interface/TrackBase.h"
00028 
00029 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00030 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
00031 
00032 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00033 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00034 #include "CircleFromThreePoints.h"
00035 
00036 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00037 #include "FWCore/Utilities/interface/InputTag.h"
00038 #include "DataFormats/GeometrySurface/interface/OpenBounds.h"
00039 
00040 /*
00041 #include "RecoVertex/KalmanVertexFit/interface/SingleTrackVertexConstraint.h"
00042 #include "TrackingTools/TransientTrack/interface/TransientTrackFromFTSFactory.h"
00043 
00044 #include "TrackingTools/TransientTrack/interface/TransientTrackFromFTSFactory.h"
00045 #include "Alignment/ReferenceTrajectories/interface/BeamSpotTransientTrackingRecHit.h"
00046 #include "Alignment/ReferenceTrajectories/interface/BeamSpotGeomDet.h"
00047 #include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
00048 #include <TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h>
00049 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToPoint.h"
00050 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToBeamLine.h"
00051 */
00052 
00053 
00054 #include <sstream>
00055 template <class T> inline T sqr( T t) {return t*t;}
00056 
00057 KFBasedPixelFitter::MyBeamSpotHit::MyBeamSpotHit (const reco::BeamSpot &beamSpot, const GeomDet * geom)
00058   : TransientTrackingRecHit(geom, DetId(0))
00059 {
00060   localPosition_ = LocalPoint(0.,0.,0.);
00061   localError_ = LocalError( sqr(beamSpot.BeamWidthX()), 0.0, sqr(beamSpot.sigmaZ())); //neglect XY differences and BS slope
00062 }
00063 
00064 AlgebraicVector KFBasedPixelFitter::MyBeamSpotHit::parameters() const 
00065 {
00066     AlgebraicVector result(1);
00067     result[0] = localPosition().x();
00068     return result;
00069 }
00070 AlgebraicSymMatrix KFBasedPixelFitter::MyBeamSpotHit::parametersError() const
00071 {
00072   LocalError le = localPositionError();
00073   AlgebraicSymMatrix m(1);
00074   m[0][0] = le.xx();
00075   return m;
00076 }
00077 AlgebraicMatrix KFBasedPixelFitter::MyBeamSpotHit::projectionMatrix() const 
00078 {
00079   AlgebraicMatrix matrix( 1, 5, 0);
00080   matrix[0][3] = 1;
00081   return matrix;
00082 }
00083 
00084 
00085 KFBasedPixelFitter::KFBasedPixelFitter( const edm::ParameterSet& cfg)
00086  :  
00087     thePropagatorLabel(cfg.getParameter<std::string>("propagator")),
00088     thePropagatorOppositeLabel(cfg.getParameter<std::string>("propagatorOpposite")),
00089     theUseBeamSpot(cfg.getParameter<bool>("useBeamSpotConstraint")),
00090     theTTRHBuilderName(cfg.getParameter<std::string>("TTRHBuilder")) 
00091 { 
00092   if (theUseBeamSpot) theBeamSpot = cfg.getParameter<edm::InputTag>("beamSpotConstraint");
00093 }
00094 
00095 reco::Track* KFBasedPixelFitter::run(
00096     const edm::Event& ev,
00097     const edm::EventSetup& es,
00098     const std::vector<const TrackingRecHit * > & hits,
00099     const TrackingRegion & region) const
00100 {
00101   int nhits = hits.size();
00102   if (nhits <2) return 0;
00103 
00104   edm::ESHandle<TrackerGeometry> tracker;
00105   es.get<TrackerDigiGeometryRecord>().get(tracker);
00106 
00107   edm::ESHandle<MagneticField> field;
00108   es.get<IdealMagneticFieldRecord>().get(field);
00109 
00110   edm::ESHandle<TransientTrackingRecHitBuilder> ttrhb;
00111   es.get<TransientRecHitRecord>().get( theTTRHBuilderName, ttrhb);
00112 
00113   float ptMin = region.ptMin();
00114 
00115   const GlobalPoint vertexPos = region.origin();
00116   GlobalError vertexErr( sqr(region.originRBound()), 0, sqr(region.originRBound()), 0, 0, sqr(region.originZBound()));
00117 
00118   std::vector<GlobalPoint> points(nhits);
00119   points[0] = tracker->idToDet(hits[0]->geographicalId())->toGlobal(hits[0]->localPosition());
00120   points[1] = tracker->idToDet(hits[1]->geographicalId())->toGlobal(hits[1]->localPosition());
00121   points[2] = tracker->idToDet(hits[2]->geographicalId())->toGlobal(hits[2]->localPosition());
00122 
00123   //
00124   //initial Kinematics
00125   //
00126   GlobalVector initMom;
00127   int charge;
00128   float theta;
00129   CircleFromThreePoints circle(points[0], points[1], points[2]);
00130   if (circle.curvature() > 1.e-4) {
00131     float invPt = PixelRecoUtilities::inversePt( circle.curvature(), es);
00132     float valPt = 1.f/invPt;
00133     float chargeTmp =    (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
00134                        - (points[1].y()-points[0].y())*(points[2].x()-points[1].x()); 
00135     int charge =  (chargeTmp>0) ? -1 : 1;
00136     float valPhi = (charge>0) ? std::atan2(circle.center().x(),-circle.center().y()) :  std::atan2(-circle.center().x(),circle.center().y());
00137     theta = GlobalVector(points[1]-points[0]).theta();
00138     initMom = GlobalVector(valPt*cos(valPhi), valPt*sin(valPhi), valPt/tan(theta)); 
00139   } 
00140   else {
00141     initMom = GlobalVector(points[1]-points[0]);
00142     initMom *= 10000./initMom.perp();
00143     charge = 1;
00144     theta = initMom.theta();
00145   }
00146   GlobalTrajectoryParameters initialKine(vertexPos, initMom, TrackCharge(charge), &*field);
00147 
00148   //
00149   // initial error
00150   //
00151   AlgebraicSymMatrix55 C = ROOT::Math::SMatrixIdentity();
00152   float sin2th = sqr(sin(theta));
00153   float minC00 = 1.0;
00154   C[0][0] = std::max(sin2th/sqr(ptMin), minC00);
00155   float zErr = vertexErr.czz();
00156   float transverseErr = vertexErr.cxx(); // assume equal cxx cyy
00157   C[3][3] = transverseErr;
00158   C[4][4] = zErr*sin2th + transverseErr*(1-sin2th);
00159   CurvilinearTrajectoryError initialError(C);
00160 
00161   FreeTrajectoryState fts(initialKine, initialError);
00162 
00163   // get propagator
00164   edm::ESHandle<Propagator>  propagator;
00165   es.get<TrackingComponentsRecord>().get(thePropagatorLabel, propagator);
00166 
00167   // get updator
00168   KFUpdator  updator;
00169 
00170   // Now update initial state track using information from hits.
00171   TrajectoryStateOnSurface outerState;
00172   DetId outerDetId = 0;
00173   const TrackingRecHit* hit = 0;
00174   for ( unsigned int iHit = 0; iHit < hits.size(); iHit++) {
00175     hit = hits[iHit];
00176     if (iHit==0) outerState = propagator->propagate(fts,tracker->idToDet(hit->geographicalId())->surface());
00177     outerDetId = hit->geographicalId();
00178     TrajectoryStateOnSurface state = propagator->propagate(outerState, tracker->idToDet(outerDetId)->surface());
00179     if (!state.isValid()) return 0;
00180 //    TransientTrackingRecHit::RecHitPointer recHit = (ttrhb->build(hit))->clone(state);
00181     TransientTrackingRecHit::RecHitPointer recHit =  ttrhb->build(hit);
00182     outerState =  updator.update(state, *recHit);
00183     if (!outerState.isValid()) return 0;
00184   }
00185   
00186 
00187 
00188   // get propagator
00189   edm::ESHandle<Propagator>  opropagator;
00190   es.get<TrackingComponentsRecord>().get(thePropagatorOppositeLabel, opropagator);
00191   TrajectoryStateOnSurface innerState = outerState;
00192   DetId innerDetId = 0;
00193   innerState.rescaleError(100000.);
00194   for ( int iHit = 2; iHit >= 0; --iHit) {
00195     hit = hits[iHit];
00196     innerDetId = hit->geographicalId();
00197     TrajectoryStateOnSurface state = opropagator->propagate(innerState, tracker->idToDet(innerDetId)->surface());
00198     if (!state.isValid()) return 0;
00199 //  TransientTrackingRecHit::RecHitPointer recHit = (ttrhb->build(hit))->clone(state);
00200     TransientTrackingRecHit::RecHitPointer recHit = ttrhb->build(hit);
00201     innerState =  updator.update(state, *recHit);
00202     if (!innerState.isValid()) return 0;
00203   }
00204 
00205 
00206   // extrapolate to vertex
00207   TrajectoryStateOnSurface  impactPointState =  TransverseImpactPointExtrapolator(&*field).extrapolate( innerState, vertexPos);
00208   if (!impactPointState.isValid()) return 0;
00209 
00210   //
00211   // optionally update impact point state with Bs constraint
00212   // using this potion makes sense if vertexPos (from TrackingRegion is centerewd at BeamSpot).
00213   //
00214   if (theUseBeamSpot) {
00215     edm::Handle<reco::BeamSpot> beamSpot;
00216     ev.getByLabel( "offlineBeamSpot", beamSpot);
00217     MyBeamSpotGeomDet bsgd(BoundPlane::build(impactPointState.surface().position(), impactPointState.surface().rotation(),OpenBounds()));
00218     MyBeamSpotHit     bsrh(*beamSpot, &bsgd);
00219     impactPointState = updator.update(impactPointState, bsrh); //update
00220     impactPointState = TransverseImpactPointExtrapolator(&*field).extrapolate( impactPointState, vertexPos); //reextrapolate
00221     if (!impactPointState.isValid()) return 0;
00222   }
00223 
00224   int ndof = 2*hits.size()-5;
00225   GlobalPoint vv = impactPointState.globalPosition();
00226   math::XYZPoint  pos( vv.x(), vv.y(), vv.z() );
00227   GlobalVector pp = impactPointState.globalMomentum();
00228   math::XYZVector mom( pp.x(), pp.y(), pp.z() );
00229 
00230   float chi2 = 0.;
00231   reco::Track * track = new reco::Track( chi2, ndof, pos, mom,
00232         impactPointState.charge(), impactPointState.curvilinearError());
00233 
00234 /*
00235     vv = outerState.globalPosition(); 
00236     pp = outerState.globalMomentum();
00237     math::XYZPoint  outerPosition( vv.x(), vv.y(), vv.z()); 
00238     math::XYZVector outerMomentum( pp.x(), pp.y(), pp.z());
00239     vv = innerState.globalPosition(); 
00240     pp = innerState.globalMomentum();
00241     math::XYZPoint  innerPosition( vv.x(), vv.y(), vv.z()); 
00242     math::XYZVector innerMomentum( pp.x(), pp.y(), pp.z());
00243 
00244     reco::TrackExtra extra( outerPosition, outerMomentum, true,
00245                       innerPosition, innerMomentum, true,
00246                       outerState.curvilinearError(), outerDetId,
00247                       innerState.curvilinearError(), innerDetId,  
00248                       anyDirection);
00249 */
00250 
00251 //  std::cout <<"TRACK CREATED" << std::endl;
00252   return track;
00253 }