CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KFFittingSmoother.h
Go to the documentation of this file.
1 
12 
15 
20 
21 
22 namespace {
23 
24  struct KFFittingSmootherParam {
25 
26  explicit KFFittingSmootherParam(const edm::ParameterSet & conf) :
27  theEstimateCut(conf.getParameter<double>("EstimateCut"))
28  ,theMaxFractionOutliers(conf.getParameter<double>("MaxFractionOutliers"))
29  ,theMaxNumberOfOutliers(conf.getParameter<int>("MaxNumberOfOutliers"))
30  ,theNoOutliersBeginEnd(conf.getParameter<bool>("NoOutliersBeginEnd"))
31  ,theMinDof(conf.getParameter<int>("MinDof"))
32  ,theMinNumberOfHits(conf.getParameter<int>("MinNumberOfHits"))
33  ,rejectTracksFlag(conf.getParameter<bool>("RejectTracks"))
34  ,breakTrajWith2ConsecutiveMissing(conf.getParameter<bool>("BreakTrajWith2ConsecutiveMissing"))
35  ,noInvalidHitsBeginEnd(conf.getParameter<bool>("NoInvalidHitsBeginEnd"))
36  {}
37 
38  double theEstimateCut;
39 
40  float theMaxFractionOutliers;
41  int theMaxNumberOfOutliers;
42  bool theNoOutliersBeginEnd;
43  int theMinDof;
44 
45  int theMinNumberOfHits;
46  bool rejectTracksFlag;
47  bool breakTrajWith2ConsecutiveMissing;
48  bool noInvalidHitsBeginEnd;
49 
50 
51  };
52 
53  class KFFittingSmoother final : public TrajectoryFitter, private KFFittingSmootherParam {
54 
55  public:
56  ~KFFittingSmoother() {}
57 
58  private:
59  KFFittingSmoother(const TrajectoryFitter& aFitter,
60  const TrajectorySmoother& aSmoother,
61  const edm::ParameterSet & conf) :
62  KFFittingSmootherParam(conf),
63  theFitter(aFitter.clone()),
64  theSmoother(aSmoother.clone())
65  {}
66 
67 
68 
69  static void fillDescriptions(edm::ParameterSetDescription & desc) {
70  desc.add<double>("EstimateCut",-1);
71  desc.add<double>("MaxFractionOutliers",0.3);
72  desc.add<int>("MaxNumberOfOutliers",3);
73  desc.add<int>("MinDof",2);
74  desc.add<bool>("NoOutliersBeginEnd",false);
75  desc.add<int>("MinNumberOfHits",5);
76  desc.add<bool>("RejectTracks",true);
77  desc.add<bool>("BreakTrajWith2ConsecutiveMissing",true);
78  desc.add<bool>("NoInvalidHitsBeginEnd",true);
79  desc.add<double>("LogPixelProbabilityCut",0);
80 
81  }
82 
83 
84  Trajectory fitOne(const Trajectory& t, fitType type) const override;
85  Trajectory fitOne(const TrajectorySeed& aSeed,
86  const RecHitContainer& hits,
87  const TrajectoryStateOnSurface& firstPredTsos, fitType type) const override;
88  Trajectory fitOne(const TrajectorySeed& aSeed,
89  const RecHitContainer& hits, fitType type) const override;
90 
91  const TrajectoryFitter* fitter() const {return theFitter.get();}
92  const TrajectorySmoother* smoother() const {return theSmoother.get();}
93 
94  std::unique_ptr<TrajectoryFitter> clone() const override {
95  return std::unique_ptr<TrajectoryFitter>(
96  new KFFittingSmoother(*theFitter,
97  *theSmoother,
98  *this)
99  );
100  }
101 
102 
103  virtual void setHitCloner(TkCloner const * hc) override {
104  theFitter->setHitCloner(hc);
105  theSmoother->setHitCloner(hc);
106  }
107 
108  KFFittingSmoother(const TrajectoryFitter& aFitter,
109  const TrajectorySmoother& aSmoother,
110  KFFittingSmootherParam const & other) :
111  KFFittingSmootherParam(other),
112  theFitter(aFitter.clone()),
113  theSmoother(aSmoother.clone())
114  {}
115 
116 
117 
118  Trajectory smoothingStep(Trajectory && fitted) const {
119  if (theEstimateCut>0) {
120  // remove "outlier" at the end of Traj
121  while (!fitted.empty() &&
122  ( !fitted.lastMeasurement().recHitR().isValid()
123  || ( fitted.lastMeasurement().recHitR().det()!=nullptr && fitted.lastMeasurement().estimate()>theEstimateCut)
124  )
125  )
126  fitted.pop();
127  if (fitted.foundHits() < theMinNumberOfHits ) return Trajectory();
128  }
129  return theSmoother->trajectory(fitted);
130  }
131 
132  private:
133 
134  const std::unique_ptr<TrajectoryFitter> theFitter;
135  const std::unique_ptr<TrajectorySmoother> theSmoother;
136 
137 
139  static bool checkForNans(const Trajectory &theTraj);
140 
141  friend class KFFittingSmootherESProducer;
142 
143  };
144 
145 
146  // #define VI_DEBUG
147 
148 #ifdef VI_DEBUG
149 #define DPRINT(x) std::cout << x << ": "
150 #define PRINT std::cout
151 #else
152 #define DPRINT(x) LogTrace(x)
153 #define PRINT LogTrace("")
154 #endif
155 
156 
157 using namespace std;
158 
159 Trajectory KFFittingSmoother::fitOne(const Trajectory& t, fitType type) const {
160  if (!t.isValid() ) return Trajectory();
161  return smoothingStep(theFitter->fitOne(t,type));
162 }
163 
164 bool KFFittingSmoother::checkForNans(const Trajectory & theTraj)
165 {
166  if (edm::isNotFinite(theTraj.chiSquared ())) return false;
167  auto const & vtm = theTraj.measurements();
168  for (auto const & tm : vtm) {
169  if (edm::isNotFinite(tm.estimate())) return false;
170  auto const & v = tm.updatedState ().localParameters ().vector();
171  for (int i=0;i<5;++i) if (edm::isNotFinite(v[i])) return false;
172  auto const & m = tm.updatedState ().curvilinearError ().matrix();
173  for (int i=0;i<5;++i)
174  for (int j=0;j<(i+1);++j) if (edm::isNotFinite(m(i,j))) return false;
175  }
176  return true;
177 }
178 
180  const RecHitContainer& hits,
181  const TrajectoryStateOnSurface& firstPredTsos,
182  fitType type) const
183 {
184  LogDebug("TrackFitters") << "In KFFittingSmoother::fit";
185 
186  if ( hits.empty() ) return Trajectory();
187 
188  RecHitContainer myHits = hits;
189  Trajectory tmp_first;
190 
191  //call the fitter
192  Trajectory smoothed = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
193 
194 
195  do {
196 
197 #ifdef EDM_ML_DEBUG
198  //if no outliers the fit is done only once
199  for (unsigned int j=0;j<myHits.size();j++) {
200  if (myHits[j]->det())
201  LogTrace("TrackFitters") << "hit #:" << j+1 << " rawId=" << myHits[j]->det()->geographicalId().rawId()
202  << " validity=" << myHits[j]->isValid();
203  else
204  LogTrace("TrackFitters") << "hit #:" << j+1 << " Hit with no Det information";
205  }
206 #endif
207 
208 
209  bool hasNaN = false;
210  if ( !smoothed.isValid() || (hasNaN = !checkForNans(smoothed)) || ( smoothed.foundHits() < theMinNumberOfHits ) ) {
211  if(hasNaN) edm::LogWarning("TrackNaN")<<"Track has NaN";
212  if ( smoothed.foundHits() < theMinNumberOfHits ) LogTrace("TrackFitters") << "smoothed.foundHits()<theMinNumberOfHits";
213  DPRINT("TrackFitters") << "smoothed invalid => trajectory rejected with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
214  if (rejectTracksFlag) {
215  return Trajectory();
216  } else {
217  std::swap(smoothed, tmp_first); // if first attempt, tmp_first would be invalid anyway
218  DPRINT("TrackFitters") << "smoothed invalid => returning orignal trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
219  }
220  break;
221  }
222 #ifdef EDM_ML_DEBUG
223  else {
224  LogTrace("TrackFitters") << "dump hits after smoothing";
225  Trajectory::DataContainer meas = smoothed[0].measurements();
226  for (Trajectory::DataContainer::iterator it=meas.begin();it!=meas.end();++it) {
227  LogTrace("TrackFitters") << "hit #" << meas.end()-it-1 << " validity=" << it->recHit()->isValid()
228  << " det=" << it->recHit()->geographicalId().rawId();
229  }
230  }
231 #endif
232 
233  if (myHits.size() !=smoothed.measurements().size())
234  DPRINT("TrackFitters") << "lost hits. before/after: " << myHits.size() <<'/' << smoothed.measurements().size()<< "\n";
235 
236  if ( theEstimateCut <= 0) break;
237 
238 
239  // Check if there are outliers
240 
241  auto msize = smoothed.measurements().size();
242  declareDynArray(unsigned int, msize, bad);
243  unsigned int nbad=0;
244  unsigned int ind=0;
245  unsigned int lastValid=smoothed.measurements().size();
246  for (auto const & tm : smoothed.measurements() ) {
247  if (tm.estimate()>theEstimateCut
248  && tm.recHitR().det()!=nullptr // do not consider outliers constraints and other special "hits"
249  ) bad[nbad++]=ind;
250  if (ind<lastValid && tm.recHitR().det()!=nullptr && tm.recHitR().isValid()) lastValid=ind;
251  ++ind;
252  }
253 
254  if (0==nbad) break;
255 
256  DPRINT("TrackFitters") << "size/found/outliers list " << smoothed.measurements().size() <<'/' << smoothed.foundHits() << ' ' << nbad << ": ";
257  for (auto i=0U; i<nbad; ++i) PRINT << bad[i] <<',';
258  PRINT << endl;
259 
260 
261 
262  if (
263  // (smoothed.foundHits() == theMinNumberOfHits) ||
264  int(nbad)>theMaxNumberOfOutliers ||
265  float(nbad) > theMaxFractionOutliers*float(smoothed.foundHits())
266  ) {
267  DPRINT("TrackFitters") << "smoothed low quality => trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
268  PRINT << "try to remove " << lastValid << std::endl;
269  nbad = 0; // try to short the traj... (below lastValid will be added)
270 
271  // do not perform outliers rejection if track is already low quality
272  /*
273  if ( rejectTracksFlag && (smoothed.chiSquared() > theEstimateCut*smoothed.ndof()) ) {
274  DPRINT("TrackFitters") << "smoothed low quality => trajectory rejected with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
275  return Trajectory();
276  } else {
277  DPRINT("TrackFitters") << "smoothed low quality => return original trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
278  }
279  break;
280  */
281  }
282 
283 
284  // always add last valid hit as outlier candidate
285  bad[nbad++]=lastValid;
286 
287  // if ( (smoothed.ndof()<theMinDof) | ) break;
288 
289 
290 
291  assert(smoothed.measurements().size() <= myHits.size());
292 
293  myHits.resize(smoothed.measurements().size()); // hits are only removed from the back...
294 
295  assert(smoothed.measurements().size() == myHits.size());
296 
297  declareDynArray(Trajectory,nbad,smoothedCand);
298 
299 
300  auto NHits = myHits.size();
302 
303  auto loc = nbad;
304  for (auto i=0U; i<nbad; ++i) {
305  auto j = NHits-bad[i]-1;
306  assert(myHits[j]->geographicalId() == smoothed.measurements()[bad[i]].recHitR().geographicalId());
307  auto removedHit = myHits[j];
308  myHits[j] = std::make_shared<InvalidTrackingRecHit>(*removedHit->det(), TrackingRecHit::missing);
309  smoothedCand[i] = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
310  myHits[j] = removedHit;
311  if (smoothedCand[i].isValid() && smoothedCand[i].chiSquared()< minChi2) {
312  minChi2 = smoothedCand[i].chiSquared();
313  loc=i;
314  }
315  }
316 
317  if (loc == nbad) {
318  DPRINT("TrackFitters") <<"New trajectories all invalid"<< "\n";
319  return Trajectory();
320  }
321 
322  DPRINT("TrackFitters") << "outlier removed " << bad[loc] << '/' << minChi2 << " was " << smoothed.chiSquared()<< "\n";
323 
324  if (minChi2>smoothed.chiSquared()) {
325 
326  DPRINT("TrackFitters") << "removing outlier makes chi2 worse " << minChi2 << '/' << smoothed.chiSquared() << "\nOri: ";
327  for (auto const & tm : smoothed.measurements() )
328  PRINT << tm.recHitR().geographicalId() << '/' << tm.estimate() <<' ';
329  PRINT << "\nNew: ";
330  for (auto const & tm : smoothedCand[loc].measurements() )
331  PRINT << tm.recHitR().geographicalId() << '/' << tm.estimate() <<' ';
332  PRINT << "\n";
333 
334  // return Trajectory();
335  // break;
336  }
337 
338  std::swap(smoothed,tmp_first);
339  myHits[NHits-bad[loc]-1] = std::make_shared<InvalidTrackingRecHit>(*myHits[NHits-bad[loc]-1]->det(), TrackingRecHit::missing);
340  std::swap(smoothed,smoothedCand[loc]);
341  // firstTry=false;
342 
343  DPRINT("TrackFitters") << "new trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
344 
345 
346 
347 
348  // Look if there are two consecutive invalid hits FIXME: take into account split matched hits!!!
349  if ( breakTrajWith2ConsecutiveMissing ) {
350  unsigned int firstinvalid = myHits.size();
351  for ( unsigned int j=0; j<myHits.size()-1; ++j )
352  {
353  if ( ((myHits[j ]->type() == TrackingRecHit::missing) && (myHits[j ]->geographicalId().rawId() != 0)) &&
354  ((myHits[j+1]->type() == TrackingRecHit::missing) && (myHits[j+1]->geographicalId().rawId() != 0))
355  && ( (myHits[j ]->geographicalId().rawId()&(~3)) != (myHits[j+1]->geographicalId().rawId()&(~3) ) ) // same gluedDet
356  )
357  {
358  firstinvalid = j;
359  DPRINT("TrackFitters") << "Found two consecutive missing hits. First invalid: " << firstinvalid << "\n";
360  break;
361  }
362  }
363 
364  //reject all the hits after the last valid before two consecutive invalid (missing) hits
365  //hits are sorted in the same order as in the track candidate FIXME??????
366  if (firstinvalid != myHits.size()) {
367  myHits.erase(myHits.begin()+firstinvalid,myHits.end());
368  smoothed = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
369  DPRINT("TrackFitters") << "Trajectory shortened " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
370  }
371  }
372 
373  } // do
374  while ( true );
375 
376  if ( smoothed.isValid() ) {
377  if ( noInvalidHitsBeginEnd
378  && !smoothed.empty() //should we send a warning ?
379  ) {
380  // discard latest dummy measurements
381  if (!smoothed.empty() && !smoothed.lastMeasurement().recHitR().isValid() )
382  LogTrace("TrackFitters") << "Last measurement is invalid";
383 
384  while (!smoothed.empty() && !smoothed.lastMeasurement().recHitR().isValid() )
385  smoothed.pop();
386 
387  //remove the invalid hits at the begin of the trajectory
388  if (!smoothed.empty() && !smoothed.firstMeasurement().recHitR().isValid() ) {
389  LogTrace("TrackFitters") << "First measurement is in`valid";
390  Trajectory tmpTraj(smoothed.seed(),smoothed.direction());
391  Trajectory::DataContainer & meas = smoothed.measurements();
392  auto it = meas.begin();
393  for ( ; it!=meas.end(); ++it )
394  if ( it->recHitR().isValid() ) break;
395  tmpTraj.push(std::move(*it),smoothed.chiSquared());//push the first valid measurement and set the same global chi2
396 
397  for (auto itt=it+1; itt!=meas.end();++itt)
398  tmpTraj.push(std::move(*itt),0);//add all the other measurements
399 
400  std::swap(smoothed,tmpTraj);
401 
402  } // if ( !smoothed[0].firstMeasurement().recHit()->isValid() )
403 
404  } // if ( noInvalidHitsBeginEnd )
405 
406  LogTrace("TrackFitters") << "end: returning smoothed trajectory with chi2="
407  << smoothed.chiSquared() ;
408 
409  //LogTrace("TrackFitters") << "dump hits before return";
410  //Trajectory::DataContainer meas = smoothed[0].measurements();
411  //for (Trajectory::DataContainer::iterator it=meas.begin();it!=meas.end();++it) {
412  //LogTrace("TrackFitters") << "hit #" << meas.end()-it-1 << " validity=" << it->recHit()->isValid()
413  //<< " det=" << it->recHit()->geographicalId().rawId();
414  //}
415 
416  }
417 
418  return smoothed;
419 
420 }
421 
422 
424  const RecHitContainer& hits,fitType type) const{
425 
426  throw cms::Exception("TrackFitters",
427  "KFFittingSmoother::fit(TrajectorySeed, <TransientTrackingRecHit>) not implemented");
428 
429  return Trajectory();
430 }
431 
432 } // namespace
#define LogDebug(id)
bool empty() const
True if trajectory has no measurements.
Definition: Trajectory.h:299
type
Definition: HCALResponse.h:21
int foundHits() const
Definition: Trajectory.h:279
int i
Definition: DBlmapReader.cc:9
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:330
tuple NHits
Definition: listHistos.py:105
assert(m_qm.get())
PropagationDirection const & direction() const
Definition: Trajectory.cc:125
DataContainer const & measurements() const
Definition: Trajectory.h:250
std::vector< TrajectoryMeasurement > DataContainer
Definition: Trajectory.h:44
bool isNotFinite(T x)
Definition: isFinite.h:10
virtual void setHitCloner(TkCloner const *)=0
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:228
def move
Definition: eostools.py:510
int j
Definition: DBlmapReader.cc:9
virtual Trajectory fitOne(const Trajectory &traj, fitType type=standard) const =0
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
bool isValid() const
Definition: Trajectory.h:324
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:241
void pop()
Definition: Trajectory.cc:16
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
ConstRecHitPointer::element_type const & recHitR() const
#define DPRINT(x)
float chiSquared() const
Definition: Trajectory.h:307
#define PRINT
susybsm::HSCParticleCollection hc
Definition: classes.h:25
#define declareDynArray(T, n, x)
Definition: DynArray.h:58