CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SingleTrackVertexConstraint.cc
Go to the documentation of this file.
5 
6 #include <algorithm>
7 using namespace std;
8 using namespace reco;
9 
10 namespace {
11  // FIXME
12  // hard-coded tracker bounds
13  // workaround while waiting for Geometry service
14  static const float TrackerBoundsRadius = 112;
15  static const float TrackerBoundsHalfLength = 273.5;
16  bool insideTrackerBounds(const GlobalPoint& point) {
17  return ((point.transverse() < TrackerBoundsRadius)
18  && (abs(point.z()) < TrackerBoundsHalfLength));
19  }
20 }
21 
23  const TransientTrack & track, const GlobalPoint& priorPos,
24  const GlobalError & priorError) const
25 {
26  VertexState priorVertexState(priorPos, priorError);
27  return constrain(track, priorVertexState);
28 }
29 
30 
32  const TransientTrack & track, const VertexState priorVertexState) const
33 {
34  // Linearize tracks
35 
36  typedef CachingVertex<5>::RefCountedVertexTrack RefCountedVertexTrack;
37  typedef VertexTrack<5>::RefCountedLinearizedTrackState RefCountedLinearizedTrackState;
38 
39  double field = track.field()->inInverseGeV(track.impactPointState().globalPosition()).z();
40  if (fabs(field) < 1e-4) {
41  LogDebug("RecoVertex/SingleTrackVertexConstraint")
42  << "Initial state is very far, field is close to zero (<1e-4): " << field << "\n";
43  return BTFtuple(false, TransientTrack(), 0.);
44  }
45 
46  RefCountedLinearizedTrackState lTrData
47  = theLTrackFactory.linearizedTrackState(priorVertexState.position(), track);
48  RefCountedVertexTrack vertexTrack = theVTrackFactory.vertexTrack(lTrData, priorVertexState);
49 
50  // Fit vertex
51 
52  std::vector<RefCountedVertexTrack> initialTracks;
53  CachingVertex<5> vertex(priorVertexState,priorVertexState,initialTracks,0);
54  vertex = vertexUpdator.add(vertex, vertexTrack);
55  if (!vertex.isValid()) {
56  return BTFtuple(false, TransientTrack(), 0.);
57  } else if (doTrackerBoundCheck_ && (!insideTrackerBounds(vertex.position()))) {
58  LogDebug("RecoVertex/SingleTrackVertexConstraint")
59  << "Fitted position is out of tracker bounds.\n";
60  return BTFtuple(false, TransientTrack(), 0.);
61  }
62 
63  RefCountedVertexTrack nTrack = theVertexTrackUpdator.update(vertex, vertexTrack);
64  return BTFtuple(true, nTrack->refittedState()->transientTrack(), nTrack->smoothedChi2());
65 }
66 
68  const FreeTrajectoryState & fts, const GlobalPoint& priorPos,
69  const GlobalError& priorError) const
70 {
71  return constrain(ttFactory.build(fts), priorPos, priorError);
72 }
73 
75  const TransientTrack & track, const reco::BeamSpot & spot ) const
76 {
77  VertexState priorVertexState(spot);
78  return constrain(track, priorVertexState);
79 }
80 
82  const FreeTrajectoryState & fts, const reco::BeamSpot & spot) const
83 {
84  VertexState priorVertexState(spot);
85  return constrain(ttFactory.build(fts), priorVertexState);
86 }
87 
#define LogDebug(id)
GlobalPoint globalPosition() const
const MagneticField * field() const
float float float z
T transverse() const
Definition: PV3DBase.h:73
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:39
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
GlobalPoint position() const
BTFtuple constrain(const reco::TransientTrack &track, const GlobalPoint &priorPos, const GlobalError &priorError) const
bool isValid() const
Definition: CachingVertex.h:96
boost::tuple< bool, reco::TransientTrack, float > BTFtuple
TrajectoryStateOnSurface impactPointState() const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5