CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SequentialVertexFitter.cc
Go to the documentation of this file.
4 
5 #include <algorithm>
6 using namespace std;
7 using namespace reco;
8 
9 namespace {
10  // FIXME
11  // hard-coded tracker bounds
12  // workaround while waiting for Geometry service
13  static const float TrackerBoundsRadius = 112;
14  static const float TrackerBoundsHalfLength = 273.5;
15  bool insideTrackerBounds(const GlobalPoint& point) {
16  return ((point.transverse() < TrackerBoundsRadius)
17  && (abs(point.z()) < TrackerBoundsHalfLength));
18  }
19 }
20 
21 
22 template <unsigned int N>
24  const LinearizationPointFinder & linP,
25  const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
26  const AbstractLTSFactory<N> & ltsf ) :
27  theLinP(linP.clone()), theUpdator(updator.clone()),
28  theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
29 {
31 }
32 
33 template <unsigned int N>
35  const edm::ParameterSet& pSet, const LinearizationPointFinder & linP,
36  const VertexUpdator<N> & updator, const VertexSmoother<N> & smoother,
37  const AbstractLTSFactory<N> & ltsf) :
38  thePSet(pSet), theLinP(linP.clone()), theUpdator(updator.clone()),
39  theSmoother(smoother.clone()), theLTrackFactory ( ltsf.clone() )
40 {
42 }
43 
44 
45 template <unsigned int N>
47  const SequentialVertexFitter & original)
48 {
49  thePSet = original.parameterSet();
50  theLinP = original.linearizationPointFinder()->clone();
51  theUpdator = original.vertexUpdator()->clone();
52  theSmoother = original.vertexSmoother()->clone();
53  theMaxShift = original.maxShift();
54  theMaxStep = original.maxStep();
55  theLTrackFactory = original.linearizedTrackStateFactory()->clone();
56 }
57 
58 
59 template <unsigned int N>
61 {
62  delete theLinP;
63  delete theUpdator;
64  delete theSmoother;
65  delete theLTrackFactory;
66 }
67 
68 
69 template <unsigned int N>
71 {
72  theMaxShift = thePSet.getParameter<double>("maxDistance"); //0.01
73  theMaxStep = thePSet.getParameter<int>("maxNbrOfIterations"); //10
74 }
75 
76 template <unsigned int N>
78 {
79  thePSet.addParameter<double>("maxDistance", 0.01);
80  thePSet.addParameter<int>("maxNbrOfIterations", 10); //10
81  readParameters();
82 }
83 
84 template <unsigned int N>
86 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks) const
87 {
88  // Linearization Point
89  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
90  if (!insideTrackerBounds(linP)) linP = GlobalPoint(0,0,0);
91 
92  // Initial vertex state, with a very large error matrix
93  AlgebraicSymMatrix we(3,1);
94  GlobalError error(we*10000);
95  VertexState state(linP, error);
96  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
97  return fit(vtContainer, state, false);
98 }
99 
100 template <unsigned int N>
102  const std::vector<RefCountedVertexTrack> & tracks,
103  const reco::BeamSpot & spot ) const
104 {
105  VertexState state(spot);
106  return fit(tracks, state, true );
107 }
108 
109 template <unsigned int N>
111 SequentialVertexFitter<N>::vertex(const std::vector<RefCountedVertexTrack> & tracks) const
112 {
113  // Initial vertex state, with a very small weight matrix
114  GlobalPoint linP = tracks[0]->linearizedTrack()->linearizationPoint();
115  AlgebraicSymMatrix we(3,1);
116  GlobalError error(we*10000);
117  VertexState state(linP, error);
118  return fit(tracks, state, false);
119 }
120 
121 
122 // Fit vertex out of a set of RecTracks.
123 // Uses the specified linearization point.
124 //
125 template <unsigned int N>
127 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks,
128  const GlobalPoint& linPoint) const
129 {
130  // Initial vertex state, with a very large error matrix
131  AlgebraicSymMatrix we(3,1);
132  GlobalError error(we*10000);
133  VertexState state(linPoint, error);
134  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
135  return fit(vtContainer, state, false);
136 }
137 
138 
143 template <unsigned int N>
145 SequentialVertexFitter<N>::vertex(const std::vector<reco::TransientTrack> & tracks,
146  const BeamSpot& beamSpot) const
147 {
148  VertexState beamSpotState(beamSpot);
149  std::vector<RefCountedVertexTrack> vtContainer;
150 
151  if (tracks.size() > 1) {
152  // Linearization Point search if there are more than 1 track
153  GlobalPoint linP = theLinP->getLinearizationPoint(tracks);
154  if (!insideTrackerBounds(linP)) linP = GlobalPoint(0,0,0);
155  AlgebraicSymMatrix we(3,1);
156  GlobalError error(we*10000);
157  VertexState lpState(linP, error);
158  vtContainer = linearizeTracks(tracks, lpState);
159  } else {
160  // otherwise take the beamspot position.
161  vtContainer = linearizeTracks(tracks, beamSpotState);
162  }
163 
164  return fit(vtContainer, beamSpotState, true);
165 }
166 
167 
168 // Fit vertex out of a set of RecTracks.
169 // Uses the position as both the linearization point AND as prior
170 // estimate of the vertex position. The error is used for the
171 // weight of the prior estimate.
172 //
173 template <unsigned int N>
175  const std::vector<reco::TransientTrack> & tracks,
176  const GlobalPoint& priorPos,
177  const GlobalError& priorError) const
178 {
179  VertexState state(priorPos, priorError);
180  std::vector<RefCountedVertexTrack> vtContainer = linearizeTracks(tracks, state);
181  return fit(vtContainer, state, true);
182 }
183 
184 // Fit vertex out of a set of VertexTracks
185 // Uses the position and error for the prior estimate of the vertex.
186 // This position is not used to relinearize the tracks.
187 //
188 template <unsigned int N>
190  const std::vector<RefCountedVertexTrack> & tracks,
191  const GlobalPoint& priorPos,
192  const GlobalError& priorError) const
193 {
194  VertexState state(priorPos, priorError);
195  return fit(tracks, state, true);
196 }
197 
198 
199 // Construct a container of VertexTrack from a set of RecTracks.
200 //
201 template <unsigned int N>
202 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack>
204  const std::vector<reco::TransientTrack> & tracks,
205  const VertexState state) const
206 {
207  GlobalPoint linP = state.position();
208  std::vector<RefCountedVertexTrack> finalTracks;
209  finalTracks.reserve(tracks.size());
210  for(vector<reco::TransientTrack>::const_iterator i = tracks.begin();
211  i != tracks.end(); i++) {
213  = theLTrackFactory->linearizedTrackState(linP, *i);
214  RefCountedVertexTrack vTrData = theVTrackFactory.vertexTrack(lTrData,state);
215  finalTracks.push_back(vTrData);
216  }
217  return finalTracks;
218 }
219 
220 
221 // Construct new a container of VertexTrack with a new linearization point
222 // and vertex state, from an existing set of VertexTrack, from which only the
223 // recTracks will be used.
224 //
225 template <unsigned int N>
226 vector<typename SequentialVertexFitter<N>::RefCountedVertexTrack>
228  const std::vector<RefCountedVertexTrack> & tracks,
229  const VertexState state) const
230 {
231 
232  GlobalPoint linP = state.position();
233  std::vector<RefCountedVertexTrack> finalTracks;
234  finalTracks.reserve(tracks.size());
235  for(typename std::vector<RefCountedVertexTrack>::const_iterator i = tracks.begin();
236  i != tracks.end(); i++) {
238  (**i).linearizedTrack()->stateWithNewLinearizationPoint(linP);
239  // RefCountedLinearizedTrackState lTrData =
240  // theLTrackFactory->linearizedTrackState(linP,
241  // (**i).linearizedTrack()->track());
242  RefCountedVertexTrack vTrData =
243  theVTrackFactory.vertexTrack(lTrData,state, (**i).weight() );
244  finalTracks.push_back(vTrData);
245  }
246  return finalTracks;
247 }
248 
249 
250 // The method where the vertex fit is actually done!
251 //
252 template <unsigned int N>
254 SequentialVertexFitter<N>::fit(const std::vector<RefCountedVertexTrack> & tracks,
255  const VertexState priorVertex, bool withPrior ) const
256 {
257  std::vector<RefCountedVertexTrack> initialTracks;
258  GlobalPoint priorVertexPosition = priorVertex.position();
259  GlobalError priorVertexError = priorVertex.error();
260 
261  CachingVertex<N> returnVertex(priorVertexPosition,priorVertexError,initialTracks,0);
262  if (withPrior) {
263  returnVertex = CachingVertex<N>(priorVertexPosition,priorVertexError,
264  priorVertexPosition,priorVertexError,initialTracks,0);
265  }
266  CachingVertex<N> initialVertex = returnVertex;
267  std::vector<RefCountedVertexTrack> globalVTracks = tracks;
268 
269  // main loop through all the VTracks
270  bool validVertex = true;
271  int step = 0;
272  GlobalPoint newPosition = priorVertexPosition;
273  GlobalPoint previousPosition;
274  do {
275  CachingVertex<N> fVertex = initialVertex;
276  // make new linearized and vertex tracks for the next iteration
277  if(step != 0) globalVTracks = reLinearizeTracks(tracks,
278  returnVertex.vertexState());
279 
280  // update sequentially the vertex estimate
281  for (typename std::vector<RefCountedVertexTrack>::const_iterator i
282  = globalVTracks.begin(); i != globalVTracks.end(); i++) {
283  fVertex = theUpdator->add(fVertex,*i);
284  if (!fVertex.isValid()) break;
285  }
286 
287  validVertex = fVertex.isValid();
288  // check tracker bounds and NaN in position
289  if (validVertex && hasNan(fVertex.position())) {
290  LogDebug("RecoVertex/SequentialVertexFitter")
291  << "Fitted position is NaN.\n";
292  validVertex = false;
293  }
294 
295  if (validVertex && !insideTrackerBounds(fVertex.position())) {
296  LogDebug("RecoVertex/SequentialVertexFitter")
297  << "Fitted position is out of tracker bounds.\n";
298  validVertex = false;
299  }
300 
301  if (!validVertex) {
302  // reset initial vertex position to (0,0,0) and force new iteration
303  // if number of steps not exceeded
305  fVertex = CachingVertex<N>(GlobalPoint(0,0,0), error,
306  initialTracks, 0);
307  }
308 
309  previousPosition = newPosition;
310  newPosition = fVertex.position();
311 
312  returnVertex = fVertex;
313  globalVTracks.clear();
314  step++;
315  } while ( (step != theMaxStep) &&
316  (((previousPosition - newPosition).transverse() > theMaxShift) ||
317  (!validVertex) ) );
318 
319  if (!validVertex) {
320  LogDebug("RecoVertex/SequentialVertexFitter")
321  << "Fitted position is invalid (out of tracker bounds or has NaN). Returned vertex is invalid\n";
322  return CachingVertex<N>(); // return invalid vertex
323  }
324 
325  if (step >= theMaxStep) {
326  LogDebug("RecoVertex/SequentialVertexFitter")
327  << "The maximum number of steps has been exceeded. Returned vertex is invalid\n";
328  return CachingVertex<N>(); // return invalid vertex
329  }
330 
331  // smoothing
332  returnVertex = theSmoother->smooth(returnVertex);
333 
334  return returnVertex;
335 }
336 
337 template class SequentialVertexFitter<5>;
338 template class SequentialVertexFitter<6>;
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::vector< RefCountedVertexTrack > reLinearizeTracks(const std::vector< RefCountedVertexTrack > &tracks, const VertexState state) const
VertexState vertexState() const
Definition: CachingVertex.h:86
const VertexUpdator< N > * vertexUpdator() const
CachingVertex< N > fit(const std::vector< RefCountedVertexTrack > &tracks, const VertexState priorVertex, bool withPrior) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define abs(x)
Definition: mlp_lapack.h:159
T transverse() const
Another name for perp()
const AbstractLTSFactory< N > * linearizedTrackStateFactory() const
GlobalPoint position() const
Definition: VertexState.h:29
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const
T transverse() const
Definition: PV3DBase.h:67
const VertexSmoother< N > * vertexSmoother() const
T z() const
Definition: PV3DBase.h:58
GlobalWeight weight() const
Definition: VertexState.h:39
const edm::ParameterSet parameterSet() const
tuple tracks
Definition: testEve_cfg.py:39
T * clone(const T *tp)
Definition: Ptr.h:42
char state
Definition: procUtils.cc:75
GlobalPoint position() const
virtual LinearizationPointFinder * clone() const =0
CLHEP::HepSymMatrix AlgebraicSymMatrix
bool isValid() const
Definition: CachingVertex.h:96
GlobalError error() const
Definition: VertexState.h:34
std::vector< RefCountedVertexTrack > linearizeTracks(const std::vector< reco::TransientTrack > &tracks, const VertexState state) const
const float maxShift() 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
const LinearizationPointFinder * linearizationPointFinder() const