CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AdaptiveVertexReconstructor.cc
Go to the documentation of this file.
7 #include <algorithm>
8 
9 using namespace std;
10 
12 {
13  vector < reco::TransientTrack > trks = old.originalTracks();
14  vector < reco::TransientTrack > newtrks;
16  static const float minweight = 1.e-8; // discard all tracks with lower weight
17  for ( vector< reco::TransientTrack >::const_iterator i=trks.begin();
18  i!=trks.end() ; ++i )
19  {
20  if ( old.trackWeight ( *i ) > minweight )
21  {
22  newtrks.push_back ( *i );
23  mp[*i]=old.trackWeight ( *i );
24  }
25  }
26 
28 
29  if ( old.hasPrior() )
30  {
31  VertexState priorstate ( old.priorPosition(), old.priorError() );
32  ret=TransientVertex ( priorstate, old.vertexState(), newtrks,
33  old.totalChiSquared(), old.degreesOfFreedom() );
34  } else {
35  ret=TransientVertex ( old.vertexState(), newtrks,
36  old.totalChiSquared(), old.degreesOfFreedom() );
37  }
38  ret.weightMap ( mp ); // set weight map
39 
40  if ( old.hasRefittedTracks() )
41  {
42  // we have refitted tracks -- copy them!
43  vector < reco::TransientTrack > newrfs;
44  vector < reco::TransientTrack > oldrfs=old.refittedTracks();
45  vector < reco::TransientTrack >::const_iterator origtrkiter=trks.begin();
46  for ( vector< reco::TransientTrack >::const_iterator i=oldrfs.begin(); i!=oldrfs.end() ; ++i )
47  {
48  if ( old.trackWeight ( *origtrkiter ) > minweight )
49  {
50  newrfs.push_back ( *i );
51  }
52  origtrkiter++;
53  }
54  if ( newrfs.size() ) ret.refittedTracks ( newrfs ); // copy refitted tracks
55  }
56 
57  if ( ret.refittedTracks().size() > ret.originalTracks().size() )
58  {
59  edm::LogError("AdaptiveVertexReconstructor" )
60  << "More refitted tracks than original tracks!";
61  }
62 
63  return ret;
64 }
65 
67  const TransientVertex & newvtx,
68  set < reco::TransientTrack > & remainingtrks,
69  float w ) const
70 {
71  /*
72  * Erase tracks that are in newvtx from remainingtrks
73  * But erase only if trackweight > w
74  */
75  const vector < reco::TransientTrack > & origtrks = newvtx.originalTracks();
76  bool erased=false;
77  for ( vector< reco::TransientTrack >::const_iterator i=origtrks.begin();
78  i!=origtrks.end(); ++i )
79  {
80  double weight = newvtx.trackWeight ( *i );
81  if ( weight > w )
82  {
83  remainingtrks.erase ( *i );
84  erased=true;
85  };
86  };
87 }
88 
90  float primcut, float seccut, float min_weight, bool smoothing ) :
91  thePrimaryFitter ( 0 ), theSecondaryFitter ( 0 ),
92  theMinWeight( min_weight ), theWeightThreshold ( 0.001 )
93 {
94  setupFitters ( primcut, 256., 0.25, seccut, 256., 0.25, smoothing );
95 }
96 
98 {
99  if ( thePrimaryFitter ) delete thePrimaryFitter;
101 }
102 
104  float primT, float primr, float seccut, float secT,
105  float secr, bool smoothing )
106 {
107  VertexSmoother<5> * smoother ;
108  if ( smoothing )
109  {
110  smoother = new KalmanVertexSmoother();
111  } else {
112  smoother = new DummyVertexSmoother<5>();
113  }
114 
115  if ( thePrimaryFitter ) delete thePrimaryFitter;
117 
118  /*
119  edm::LogError ("AdaptiveVertexReconstructor" )
120  << "Tini and r are hardcoded now!";
121  */
125  // if the primary fails, sth is wrong, so here we set a threshold on the weight.
130  // need to set it or else we have
131  // unwanted exceptions to deal with.
132  // cleanup can come later!
133  delete smoother;
134 }
135 
137  : thePrimaryFitter(0), theSecondaryFitter(0), theMinWeight(0.5), theWeightThreshold ( 0.001 )
138 {
139  float primcut = 2.0;
140  float seccut = 6.0;
141  bool smoothing=false;
142  // float primT = 4096.;
143  // float primr = 0.125;
144  float primT = 256.;
145  float primr = 0.25;
146  float secT = 256.;
147  float secr = 0.25;
148 
149  try {
150  primcut = m.getParameter<double>("primcut");
151  primT = m.getParameter<double>("primT");
152  primr = m.getParameter<double>("primr");
153  seccut = m.getParameter<double>("seccut");
154  secT = m.getParameter<double>("secT");
155  secr = m.getParameter<double>("secr");
156  theMinWeight = m.getParameter<double>("minweight");
157  theWeightThreshold = m.getParameter<double>("weightthreshold");
158  smoothing = m.getParameter<bool>("smoothing");
159  } catch ( edm::Exception & e ) {
160  edm::LogError ("AdaptiveVertexReconstructor") << e.what();
161  }
162 
163  setupFitters ( primcut, primT, primr, seccut, secT, secr, smoothing );
164 }
165 
166 vector<TransientVertex>
167  AdaptiveVertexReconstructor::vertices(const vector<reco::TransientTrack> & t,
168  const reco::BeamSpot & s ) const
169 {
170  return vertices ( vector<reco::TransientTrack>(), t, s, false, true );
171 }
172 
173 vector<TransientVertex>
174  AdaptiveVertexReconstructor::vertices(const vector<reco::TransientTrack> & primaries,
175  const vector<reco::TransientTrack > & tracks, const reco::BeamSpot & s ) const
176 {
177  return vertices ( primaries, tracks, s, true, true );
178 }
179 
180 
181 vector<TransientVertex> AdaptiveVertexReconstructor::vertices (
182  const vector<reco::TransientTrack> & tracks ) const
183 {
184  return vertices ( vector<reco::TransientTrack>(), tracks, reco::BeamSpot(),
185  false, false );
186 }
187 
188 vector<TransientVertex> AdaptiveVertexReconstructor::vertices (
189  const vector<reco::TransientTrack> & primaries,
190  const vector<reco::TransientTrack> & tracks,
191  const reco::BeamSpot & s, bool has_primaries, bool usespot ) const
192 {
193  vector < TransientVertex > ret;
194  set < reco::TransientTrack > remainingtrks;
195 
196  copy(tracks.begin(), tracks.end(),
197  inserter(remainingtrks, remainingtrks.begin()));
198 
199  int ctr=0;
200  unsigned int n_tracks = remainingtrks.size();
201 
202  // cout << "[AdaptiveVertexReconstructor] DEBUG ::vertices!!" << endl;
203  try {
204  while ( remainingtrks.size() > 1 )
205  {
206  /*
207  cout << "[AdaptiveVertexReconstructor] next round: "
208  << remainingtrks.size() << endl;
209  */
210  ctr++;
211  const AdaptiveVertexFitter * fitter = theSecondaryFitter;
212  if ( ret.size() == 0 )
213  {
214  fitter = thePrimaryFitter;
215  };
216  vector < reco::TransientTrack > fittrks;
217  fittrks.reserve ( remainingtrks.size() );
218 
219  copy(remainingtrks.begin(), remainingtrks.end(), back_inserter(fittrks));
220 
221  TransientVertex tmpvtx;
222  if ( (ret.size() == 0) && has_primaries )
223  {
224  // add the primaries to the fitted tracks.
225  copy ( primaries.begin(), primaries.end(), back_inserter(fittrks) );
226  }
227  if ( (ret.size() == 0) && usespot )
228  {
229  tmpvtx=fitter->vertex ( fittrks, s );
230  } else {
231  tmpvtx=fitter->vertex ( fittrks );
232  }
233  TransientVertex newvtx = cleanUp ( tmpvtx );
234  ret.push_back ( newvtx );
235  erase ( newvtx, remainingtrks, theMinWeight );
236  if ( n_tracks == remainingtrks.size() )
237  {
238  if ( usespot )
239  {
240  // try once more without beamspot constraint!
241  usespot=false;
242  LogDebug("AdaptiveVertexReconstructor")
243  << "no tracks in vertex. trying again without beamspot constraint!";
244  continue;
245  }
246  LogDebug("AdaptiveVertexReconstructor") << "all tracks (" << n_tracks
247  << ") would be recycled for next fit. Trying with low threshold!";
248  erase ( newvtx, remainingtrks, 1.e-5 );
249  if ( n_tracks == remainingtrks.size() )
250  {
251  LogDebug("AdaptiveVertexReconstructor") << "low threshold didnt help! "
252  << "Discontinue procedure!";
253  break;
254  }
255  };
256 
257  // cout << "[AdaptiveVertexReconstructor] erased" << endl;
258  n_tracks = remainingtrks.size();
259  };
260  } catch ( VertexException & v ) {
261  // Will catch all (not enough significant tracks exceptions.
262  // in this case, the iteration can safely terminate.
263  };
264 
265  return cleanUpVertices ( ret );
266 }
267 
269  const vector < TransientVertex > & old ) const
270 {
271  vector < TransientVertex > ret;
272  for ( vector< TransientVertex >::const_iterator i=old.begin(); i!=old.end() ; ++i )
273  {
274  if (!(i->hasTrackWeight()))
275  { // if we dont have track weights, we take the vtx
276  ret.push_back ( *i );
277  continue;
278  }
279 
280  // maybe this should be replaced with asking for the ndf ...
281  // e.g. if ( ndf > - 1. )
282  int nsig=0; // number of significant tracks.
284  for ( TransientVertex::TransientTrackToFloatMap::const_iterator w=wm.begin(); w!=wm.end() ; ++w )
285  {
286  if (w->second > theWeightThreshold) nsig++;
287  }
288  if ( nsig > 1 ) ret.push_back ( *i );
289  }
290 
291  return ret;
292 }
#define LogDebug(id)
virtual char const * what() const
Definition: Exception.cc:97
AdaptiveVertexReconstructor(float primcut=2.0, float seccut=6.0, float minweight=0.5, bool smoothing=false)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
GlobalError priorError() const
TransientVertex cleanUp(const TransientVertex &old) const
Common base class.
float totalChiSquared() const
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
bool hasRefittedTracks() const
std::vector< reco::TransientTrack > originalTracks() const
void setupFitters(float primcut, float primT, float primr, float seccut, float secT, float secr, bool smoothing)
float degreesOfFreedom() const
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &v) const
bool hasPrior() const
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &) const
std::vector< TransientVertex > cleanUpVertices(const std::vector< TransientVertex > &) const
VertexState vertexState() const
AdaptiveVertexFitter * theSecondaryFitter
float trackWeight(const reco::TransientTrack &track) const
std::vector< reco::TransientTrack > refittedTracks() const
tuple tracks
Definition: testEve_cfg.py:39
GlobalPoint priorPosition() const
string s
Definition: asciidump.py:422
void erase(const TransientVertex &newvtx, std::set< reco::TransientTrack > &remainingtrks, float w) const
mathSSE::Vec4< T > v
TransientTrackToFloatMap weightMap() const