CMS 3D CMS Logo

SequentialVertexFitter.cc
Go to the documentation of this file.
4 
5 #include <algorithm>
6 using namespace std;
7 using namespace reco;
8 
9 template <unsigned int N>
12  const VertexSmoother<N>& smoother,
13  const AbstractLTSFactory<N>& ltsf)
14  : theLinP(linP.clone()),
15  theUpdator(updator.clone()),
16  theSmoother(smoother.clone()),
17  theLTrackFactory(ltsf.clone()) {
19 }
20 
21 template <unsigned int N>
23  const LinearizationPointFinder& linP,
25  const VertexSmoother<N>& smoother,
26  const AbstractLTSFactory<N>& ltsf)
27  : thePSet(pSet),
28  theLinP(linP.clone()),
29  theUpdator(updator.clone()),
30  theSmoother(smoother.clone()),
31  theLTrackFactory(ltsf.clone()) {
33 }
34 
35 template <unsigned int N>
37  thePSet = original.parameterSet();
38  theLinP = original.linearizationPointFinder()->clone();
39  theUpdator = original.vertexUpdator()->clone();
40  theSmoother = original.vertexSmoother()->clone();
41  theMaxShift = original.maxShift();
42  theMaxStep = original.maxStep();
43  theLTrackFactory = original.linearizedTrackStateFactory()->clone();
44 }
45 
46 template <unsigned int N>
48  delete theLinP;
49  delete theUpdator;
50  delete theSmoother;
51  delete theLTrackFactory;
52 }
53 
54 template <unsigned int N>
56  theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
57  theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
58 }
59 
60 template <unsigned int N>
62  thePSet.addParameter<double>("maxDistance", 0.01);
63  thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
64  readParameters();
65 }
66 
67 template <unsigned int N>
68 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack>& tracks) const {
69  // Linearization Point
70  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
71  if (!insideTrackerBounds(linP))
72  linP = GlobalPoint(0, 0, 0);
73 
74  // Initial vertex state, with a very large error matrix
75  ROOT::Math::SMatrixIdentity id;
76  AlgebraicSymMatrix33 we(id);
77  GlobalError error(we * 10000);
78  VertexState state(linP, error);
79  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
80  return fit(vtContainer, state, false);
81 }
82 
83 template <unsigned int N>
84 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack>& tracks,
85  const reco::BeamSpot& spot) const {
86  VertexState state(spot);
87  return fit(tracks, state, true);
88 }
89 
90 template <unsigned int N>
91 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack>& tracks) const {
92  // Initial vertex state, with a very small weight matrix
93  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
94  ROOT::Math::SMatrixIdentity id;
95  AlgebraicSymMatrix33 we(id);
96  GlobalError error(we * 10000);
97  VertexState state(linP, error);
98  return fit(tracks, state, false);
99 }
100 
101 // Fit vertex out of a set of RecTracks.
102 // Uses the specified linearization point.
103 //
104 template <unsigned int N>
105 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack>& tracks,
106  const GlobalPoint& linPoint) const {
107  // Initial vertex state, with a very large error matrix
108  ROOT::Math::SMatrixIdentity id;
109  AlgebraicSymMatrix33 we(id);
110  GlobalError error(we * 10000);
111  VertexState state(linPoint, error);
112  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
113  return fit(vtContainer, state, false);
114 }
115 
120 template <unsigned int N>
121 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack>& tracks,
122  const BeamSpot& beamSpot) const {
123  VertexState beamSpotState(beamSpot);
124  std::vector<RefCountedVertexTrack> vtContainer;
125 
126  if (tracks.size() > 1) {
127  // Linearization Point search if there are more than 1 track
128  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
129  if (!insideTrackerBounds(linP))
130  linP = GlobalPoint(0, 0, 0);
131  ROOT::Math::SMatrixIdentity id;
132  AlgebraicSymMatrix33 we(id);
133  GlobalError error(we * 10000);
134  VertexState lpState(linP, error);
135  vtContainer = linearizeTracks(tracks, lpState);
136  } else {
137  // otherwise take the beamspot position.
138  vtContainer = linearizeTracks(tracks, beamSpotState);
139  }
140 
141  return fit(vtContainer, beamSpotState, true);
142 }
143 
144 // Fit vertex out of a set of RecTracks.
145 // Uses the position as both the linearization point AND as prior
146 // estimate of the vertex position. The error is used for the
147 // weight of the prior estimate.
148 //
149 template <unsigned int N>
150 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack>& tracks,
151  const GlobalPoint& priorPos,
152  const GlobalError& priorError) const {
153  VertexState state(priorPos, priorError);
154  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
155  return fit(vtContainer, state, true);
156 }
157 
158 // Fit vertex out of a set of VertexTracks
159 // Uses the position and error for the prior estimate of the vertex.
160 // This position is not used to relinearize the tracks.
161 //
162 template <unsigned int N>
163 CachingVertex<N> SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack>& tracks,
164  const GlobalPoint& priorPos,
165  const GlobalError& priorError) const {
166  VertexState state(priorPos, priorError);
167  return fit(tracks, state, true);
168 }
169 
170 // Construct a container of VertexTrack from a set of RecTracks.
171 //
172 template <unsigned int N>
173 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> SequentialVertexFitter<N>::linearizeTracks(
174  const std::vector<reco::TransientTrack>& tracks, const VertexState state) const {
175  GlobalPoint linP = state.position();
176  std::vector<RefCountedVertexTrack> finalTracks;
177  finalTracks.reserve(tracks.size());
178  for (vector<reco::TransientTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
179  RefCountedLinearizedTrackState lTrData = theLTrackFactory->linearizedTrackState(linP, *i);
180  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state);
181  finalTracks.push_back(vTrData);
182  }
183  return finalTracks;
184 }
185 
186 // Construct new a container of VertexTrack with a new linearization point
187 // and vertex state, from an existing set of VertexTrack, from which only the
188 // recTracks will be used.
189 //
190 template <unsigned int N>
191 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack> SequentialVertexFitter<N>::reLinearizeTracks(
192  const std::vector<RefCountedVertexTrack>& tracks, const VertexState state) const {
193  GlobalPoint linP = state.position();
194  std::vector<RefCountedVertexTrack> finalTracks;
195  finalTracks.reserve(tracks.size());
196  for (typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin(); i != tracks.end(); i++) {
197  RefCountedLinearizedTrackState lTrData = (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
198  // RefCountedLinearizedTrackState lTrData =
199  // theLTrackFactory->linearizedTrackState(linP,
200  // (**i).linearizedTrack()->track());
201  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData, state, (**i).weight());
202  finalTracks.push_back(vTrData);
203  }
204  return finalTracks;
205 }
206 
207 // The method where the vertex fit is actually done!
208 //
209 template <unsigned int N>
210 CachingVertex<N> SequentialVertexFitter<N>::fit(const std::vector<RefCountedVertexTrack>& tracks,
211  const VertexState priorVertex,
212  bool withPrior) const {
213  std::vector<RefCountedVertexTrack> initialTracks;
214  GlobalPoint priorVertexPosition = priorVertex.position();
215  GlobalError priorVertexError = priorVertex.error();
216 
217  CachingVertex<N> returnVertex(priorVertexPosition, priorVertexError, initialTracks, 0);
218  if (withPrior) {
219  returnVertex = CachingVertex<N>(
220  priorVertexPosition, priorVertexError, priorVertexPosition, priorVertexError, initialTracks, 0);
221  }
222  CachingVertex<N> initialVertex = returnVertex;
223  std::vector<RefCountedVertexTrack> globalVTracks = tracks;
224 
225  // main loop through all the VTracks
226  bool validVertex = true;
227  int step = 0;
228  GlobalPoint newPosition = priorVertexPosition;
229  GlobalPoint previousPosition;
230  do {
231  CachingVertex<N> fVertex = initialVertex;
232  // make new linearized and vertex tracks for the next iteration
233  if (step != 0)
234  globalVTracks = reLinearizeTracks(tracks, returnVertex.vertexState());
235 
236  // update sequentially the vertex estimate
237  for (typename std::vector<RefCountedVertexTrack>::const_iterator i = globalVTracks.begin();
238  i != globalVTracks.end();
239  i++) {
240  fVertex = theUpdator->add(fVertex, *i);
241  if (!fVertex.isValid())
242  break;
243  }
244 
245  validVertex = fVertex.isValid();
246  // check tracker bounds and NaN in position
247  if (validVertex && hasNan(fVertex.position())) {
248  LogDebug("RecoVertex/SequentialVertexFitter") << "Fitted position is NaN.\n";
249  validVertex = false;
250  }
251 
252  if (validVertex && !insideTrackerBounds(fVertex.position())) {
253  LogDebug("RecoVertex/SequentialVertexFitter") << "Fitted position is out of tracker bounds.\n";
254  validVertex = false;
255  }
256 
257  if (!validVertex) {
258  // reset initial vertex position to (0,0,0) and force new iteration
259  // if number of steps not exceeded
260  ROOT::Math::SMatrixIdentity id;
261  AlgebraicSymMatrix33 we(id);
262  GlobalError error(we * 10000);
263  fVertex = CachingVertex<N>(GlobalPoint(0, 0, 0), error, initialTracks, 0);
264  }
265 
266  previousPosition = newPosition;
267  newPosition = fVertex.position();
268 
269  returnVertex = fVertex;
270  globalVTracks.clear();
271  step++;
272  } while ((step != theMaxStep) && (((previousPosition - newPosition).transverse() > theMaxShift) || (!validVertex)));
273 
274  if (!validVertex) {
275  LogDebug("RecoVertex/SequentialVertexFitter")
276  << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
277  return CachingVertex<N>(); // return invalid vertex
278  }
279 
280  if (step >= theMaxStep) {
281  LogDebug("RecoVertex/SequentialVertexFitter")
282  << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
283  return CachingVertex<N>(); // return invalid vertex
284  }
285 
286  // smoothing
287  returnVertex = theSmoother->smooth(returnVertex);
288 
289  return returnVertex;
290 }
291 
292 template class SequentialVertexFitter<5>;
293 template class SequentialVertexFitter<6>;
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
VertexState const & vertexState() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
T transverse() const
Another name for perp()
CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const override
GlobalError error() const
Definition: VertexState.h:64
bool isValid() const
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
fixed size matrix
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
step
Definition: StallMonitor.cc:98
GlobalPoint position() const
#define LogDebug(id)
GlobalPoint position() const
Definition: VertexState.h:62