CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoVertex/KalmanVertexFit/src/SingleTrackVertexConstraint.cc

Go to the documentation of this file.
00001 #include "RecoVertex/KalmanVertexFit/interface/SingleTrackVertexConstraint.h"
00002 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "MagneticField/Engine/interface/MagneticField.h" 
00005 
00006 #include <algorithm>
00007 using namespace std;
00008 using namespace reco;
00009 
00010 namespace {
00011   // FIXME
00012   // hard-coded tracker bounds
00013   // workaround while waiting for Geometry service
00014   static const float TrackerBoundsRadius = 112;
00015   static const float TrackerBoundsHalfLength = 273.5;
00016   bool insideTrackerBounds(const GlobalPoint& point) {
00017     return ((point.transverse() < TrackerBoundsRadius)
00018         && (abs(point.z()) < TrackerBoundsHalfLength));
00019   }
00020 }
00021 
00022 SingleTrackVertexConstraint::BTFtuple SingleTrackVertexConstraint::constrain(
00023         const TransientTrack & track, const GlobalPoint& priorPos,
00024         const GlobalError & priorError) const
00025 { 
00026   VertexState priorVertexState(priorPos, priorError);
00027   return constrain(track, priorVertexState);
00028 }
00029 
00030 
00031 SingleTrackVertexConstraint::BTFtuple SingleTrackVertexConstraint::constrain(
00032         const TransientTrack & track,  const VertexState priorVertexState) const
00033 {
00034   // Linearize tracks
00035 
00036   typedef CachingVertex<5>::RefCountedVertexTrack RefCountedVertexTrack;
00037   typedef VertexTrack<5>::RefCountedLinearizedTrackState RefCountedLinearizedTrackState;
00038 
00039   double field  = track.field()->inInverseGeV(track.impactPointState().globalPosition()).z();
00040   if (fabs(field) < 1e-4) {
00041       LogDebug("RecoVertex/SingleTrackVertexConstraint") 
00042          << "Initial state is very far, field is close to zero (<1e-4): " << field << "\n";
00043       return BTFtuple(false, TransientTrack(), 0.);
00044   }
00045 
00046   RefCountedLinearizedTrackState lTrData 
00047       = theLTrackFactory.linearizedTrackState(priorVertexState.position(), track);
00048   RefCountedVertexTrack vertexTrack =  theVTrackFactory.vertexTrack(lTrData, priorVertexState);
00049 
00050   // Fit vertex
00051 
00052   std::vector<RefCountedVertexTrack> initialTracks;
00053   CachingVertex<5> vertex(priorVertexState,priorVertexState,initialTracks,0);
00054   vertex = vertexUpdator.add(vertex, vertexTrack);
00055   if (!vertex.isValid()) {
00056     return BTFtuple(false, TransientTrack(), 0.);
00057   } else  if (doTrackerBoundCheck_ && (!insideTrackerBounds(vertex.position()))) {
00058       LogDebug("RecoVertex/SingleTrackVertexConstraint") 
00059          << "Fitted position is out of tracker bounds.\n";
00060       return BTFtuple(false, TransientTrack(), 0.);
00061   }
00062 
00063   RefCountedVertexTrack nTrack = theVertexTrackUpdator.update(vertex, vertexTrack);
00064   return BTFtuple(true, nTrack->refittedState()->transientTrack(), nTrack->smoothedChi2());
00065 }
00066 
00067 SingleTrackVertexConstraint::BTFtuple SingleTrackVertexConstraint::constrain(
00068         const FreeTrajectoryState & fts, const GlobalPoint& priorPos,
00069         const GlobalError& priorError) const
00070 { 
00071   return constrain(ttFactory.build(fts), priorPos, priorError);
00072 }
00073 
00074 SingleTrackVertexConstraint::BTFtuple SingleTrackVertexConstraint::constrain(
00075         const TransientTrack & track, const reco::BeamSpot & spot ) const
00076 {
00077   VertexState priorVertexState(spot);
00078   return constrain(track, priorVertexState);
00079 }
00080 
00081 SingleTrackVertexConstraint::BTFtuple SingleTrackVertexConstraint::constrain(
00082         const FreeTrajectoryState & fts, const reco::BeamSpot & spot) const
00083 { 
00084   VertexState priorVertexState(spot);
00085   return constrain(ttFactory.build(fts), priorVertexState);
00086 }
00087