StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
FwdTracker.h
1 #ifndef FWD_TRACKER_H
2 #define FWD_TRACKER_H
3 
4 #include "Fit/Fitter.h"
5 #include "TFile.h"
6 #include "TGraph2D.h"
7 #include "TH1.h"
8 #include "TH1F.h"
9 #include "TH2F.h"
10 #include "TRandom3.h"
11 #include "TTree.h"
12 #include "TVector3.h"
13 #include "TLorentzVector.h"
14 
15 #include <algorithm>
16 #include <memory>
17 #include <string>
18 #include <vector>
19 #include <unordered_map>
20 #include <fstream>
21 
22 #include "StFwdTrackMaker/include/Tracker/FwdHit.h"
23 #include "StFwdTrackMaker/include/Tracker/FwdDataSource.h"
24 #include "StFwdTrackMaker/include/Tracker/TrackFitter.h"
25 
26 #include "Criteria/Criteria.h"
27 #include "Criteria/ICriterion.h"
28 #include "KiTrack/Automaton.h"
29 #include "KiTrack/SegmentBuilder.h"
30 #include "KiTrack/SubsetHopfieldNN.h"
31 
32 #include "StFwdTrackMaker/include/Tracker/CriteriaKeeper.h"
33 
34 #include "GenFit/FitStatus.h"
35 #include "GenFit/GFRaveVertexFactory.h"
36 
37 #include "StFwdTrackMaker/FwdTrackerConfig.h"
38 #include "StFwdTrackMaker/Common.h"
39 
40 // Utility class for evaluating ID and QA truth
41 struct MCTruthUtils {
42 
43  static int dominantContribution(Seed_t hits, double &qa) {
44 
45  // track_id, hits on track
46  std::unordered_map<int,int> truth;
47  for ( auto hit : hits ) {
48  FwdHit* fhit = dynamic_cast<FwdHit*>(hit);
49  truth[ fhit->_tid ]++;
50  }
51 
52  if ( truth.size() == 0 ){
53  return -1;
54  }
55 
56  using namespace std;
57  using P = decltype(truth)::value_type;
58  auto dom = max_element(begin(truth), end(truth), [](P a, P b){ return a.second < b.second; });
59 
60  // QA represents the percentage of hits which
61  // vote the same way on the track
62  if ( hits.size() > 0 )
63  qa = double(dom->second) / double(hits.size()) ;
64  else
65  qa = 0;
66 
67  return dom->first;
68  };
69 };
70 
71 class EventStats {
72  public:
73  TString classname() const { return "EventStats"; }
74  void reset(){
75  mNumSeeds = 0;
76  mAttemptedFits = 0;
77  mFailedFits = 0;
78  mGoodFits = 0;
79  mGoodCardinals = 0;
80 
81  mAttemptedReFits = 0;
82  mFailedReFits = 0;
83  mGoodReFits = 0;
84  mPossibleReFit = 0;
85 
86  mNumFwdVertices = 0;
87  mAttemptedPrimaryFits = 0;
88  mGoodPrimaryFits = 0;
89  mFailedPrimaryFits = 0;
90 
91  mStep1Duration.clear();
92  mStep2Duration.clear();
93  mStep3Duration.clear();
94  mStep4Duration.clear();
95  mFitDuration.clear();
96  }
97  int mNumSeeds = 0;
98  int mAttemptedFits = 0;
99  int mFailedFits = 0;
100  int mGoodFits = 0;
101  int mGoodCardinals = 0;
102 
103  int mAttemptedReFits = 0;
104  int mFailedReFits = 0;
105  int mGoodReFits = 0;
106  int mPossibleReFit = 0;
107 
108  int mNumFwdVertices = 0;
109  int mAttemptedPrimaryFits = 0;
110  int mGoodPrimaryFits = 0;
111  int mFailedPrimaryFits = 0;
112  vector<float> mStep1Duration;
113  vector<float> mSeedFindingDuration;
114  vector<float> mStep2Duration;
115  vector<float> mStep3Duration;
116  vector<float> mStep4Duration;
117  vector<float> mFitDuration;
118 };
119 
121 public:
123  GenfitTrackResult( Seed_t &seed, std::shared_ptr<genfit::Track> track ) {
124  set( seed, track );
125  }
127  // Clear();
128  }
129  void Clear() {
130  if ( mTrack ){
131  mTrack->Clear();
132  }
133  }
134  void set( Seed_t &seeds, std::shared_ptr<genfit::Track> track ){
135  setSeed( seeds );
136  setTrack( track );
137  }
138  void setSeed( Seed_t &seed ){
139  mSeed = seed;
140  mIdTruth = MCTruthUtils::dominantContribution( seed, mQaTruth );
141  LOG_INFO << "GenFitTrackResult::mIdTruth = " << mIdTruth << ", QaTruth = " << mQaTruth << endm;
142  }
143  void setTrack( std::shared_ptr<genfit::Track> track ){
144  try {
145  // this->track = new genfit::Track(*track);
146  mTrack = track;
147  mTrack ->setMcTrackId(mIdTruth);
148  mStatus = mTrack->getFitStatus();
149  mTrackRep = mTrack->getCardinalRep();
150 
151  mIsFitConverged = mStatus->isFitConverged();
152  mIsFitConvergedFully = mStatus->isFitConvergedFully();
153  mIsFitConvergedPartially = mStatus->isFitConvergedPartially();
154  mNFailedPoints = mStatus->getNFailedPoints();
155  mCharge = mStatus->getCharge();
156  mChi2 = mStatus->getChi2();
157 
158  if ( mIsFitConverged ){
159  LOG_INFO << "GTR Setting momentum from track" << endm;
160  mMomentum = mTrackRep->getMom( mTrack->getFittedState(0, mTrackRep) );
161  }
162  LOG_DEBUG << "GenfitTrackResult::set Track successful" << endm;
163  } catch ( genfit::Exception &e ) {
164  LOG_ERROR << "Unable to set track -> GenfitException: " << e.what() << endm;
165  this->mTrack = nullptr;
166  this->mTrackRep = nullptr;
167 
168  this->mIsFitConverged = false;
169  this->mIsFitConvergedFully = false;
170  this->mIsFitConvergedPartially = false;
171  this->mNFailedPoints = 99;
172  this->mCharge = 0;
173  this->mChi2 = -1;
174  }
175  }
176 
180  void setDCA( TVector3 pv ){
181  mPV = pv;
182  if ( mTrack ){
183  try {
184  auto dcaState = mTrack->getFittedState( 0 );
185  this->mTrackRep->extrapolateToPoint( dcaState, mPV );
186  this->mDCA = dcaState.getPos();
187  } catch ( genfit::Exception &e ) {
188  LOG_ERROR << "CANNOT GET DCA : GenfitException: " << e.what() << endm;
189  this->mDCA = TVector3(99,99,99);
190  }
191 
192  }
193  }
194  size_t numFtt() const {
195  size_t n = 0;
196  for ( auto hit : mSeed ){
197  if ( dynamic_cast<FwdHit*>(hit)->_detid == kFttId ){
198  n++;
199  }
200  }
201  return n;
202  }
203  size_t numFst() const {
204  size_t n = 0;
205  for ( auto hit : mSeed ){
206  if ( dynamic_cast<FwdHit*>(hit)->_detid == kFstId ){
207  n++;
208  }
209  }
210  return n;
211  }
212  size_t numPV() const {
213  size_t n = 0;
214  for ( auto hit : mSeed ){
215  if ( dynamic_cast<FwdHit*>(hit)->isPV() ){
216  n++;
217  }
218  }
219  return n;
220  }
221 
222  void mergeSeeds( GenfitTrackResult &other ){
223  // combine the unique Ftt and Fst seeds
224  for ( auto hit : other.mSeed ){
225  if ( std::find( mSeed.begin(), mSeed.end(), hit ) == mSeed.end() ){
226  mSeed.push_back( hit );
227  }
228  }
229  }
230 
231  bool isPrimary = false;
232  Seed_t mSeed;
233  TVector3 mPV; // as a TVector3
234  TVector3 mMomentum;
235  float mCharge = 0;
236  float mChi2 = -1;
237  genfit::FitStatus *mStatus = nullptr;
238  genfit::AbsTrackRep *mTrackRep = nullptr;
239  std::shared_ptr<genfit::Track> mTrack = nullptr;
240  bool mIsFitConverged = false;
241  bool mIsFitConvergedFully = false;
242  bool mIsFitConvergedPartially = false;
243  size_t mNFailedPoints = 0;
244  TVector3 mDCA;
245  int mIdTruth = -1;
246  double mQaTruth = 0;
247 }; // GenfitTrackResult
248 
250  public:
251  ForwardTrackMaker() : mConfigFile("config.xml"), mEventVertex(-999, -999, -999) {
252  // noop
253  }
254 
255  const std::vector<GenfitTrackResult> &getTrackResults() const { return mTrackResults; }
256  const std::vector<Seed_t> &getTrackSeeds() const { return mTrackSeeds; }
257  const EventStats &getEventStats() const { return mEventStats; }
258 
259  void Clear(){
260  for ( auto gtr : mTrackResults ){
261  gtr.Clear();
262  }
263  mTrackResults.clear();
264  }
265 
271  void setConfigFile(std::string cf) {
272  mConfigFile = cf;
273  }
274 
280  void setSaveCriteriaValues(bool save) {
281  mSaveCriteriaValues = save;
282  }
283 
284  // Adopt external configuration file
285  void setConfig(FwdTrackerConfig cfg) { mConfig = cfg; }
286  // Adopt external hit loader
287  void setData(std::shared_ptr<FwdDataSource>data) { mDataSource = data; }
288 
295  virtual void initialize( TString geoCache, bool genHistograms) {
296  mGeoCache = geoCache;
297  mDoTrackFitting = mConfig.get<bool>("TrackFitter:active", true);
298 
299  if (!mConfig.exists("TrackFitter"))
300  mDoTrackFitting = false;
301  } //initialize
302 
309  std::vector<KiTrack::ICriterion *> loadCriteria(string path) {
310 
311  std::vector<KiTrack::ICriterion *> crits;
312  auto paths = mConfig.childrenOf(path);
313 
314  for (string p : paths) {
315  string name = mConfig.get<string>(p + ":name", "");
316  bool active = mConfig.get<bool>(p + ":active", true);
317 
318  if (false == active) {
319  continue;
320  }
321 
322  float vmin = mConfig.get<float>(p + ":min", 0);
323  float vmax = mConfig.get<float>(p + ":max", 1);
324 
325  KiTrack::ICriterion * crit = nullptr;
326  if ( name == "Crit2_BDT" ){
327  // crit = new BDTCrit2( vmin, vmax );
328  LOG_WARN << "BDT Criteria not implemented/out of date" << endm;
329  } else {
330  crit = KiTrack::Criteria::createCriterion(name, vmin, vmax);
331  }
332 
333  crit->setSaveValues(mSaveCriteriaValues);
334 
335  if (mSaveCriteriaValues)
336  crits.push_back(new CriteriaKeeper(crit)); // CriteriaKeeper intercepts values and saves them
337  else
338  crits.push_back(crit);
339 
340  }
341 
342  return crits;
343  } // loadCriteria
344 
349  void clearCriteria( std::vector<KiTrack::ICriterion *> &crits ){
350  for ( size_t i = 0; i < crits.size(); i++ ){
351  if ( crits[i] ){
352  delete crits[i];
353  crits[i] = nullptr;
354  }
355  }
356  crits.clear();
357  }
358 
365  std::vector<float> getCriteriaValues(std::string crit_name) {
366  std::vector<float> em;
367  if (mSaveCriteriaValues != true) {
368  return em;
369  }
370 
371  for (auto crit : mTwoHitCrit) {
372  if (crit_name == crit->getName()) {
373  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
374  return critKeeper->getValues();
375  }
376  }
377 
378  for (auto crit : mThreeHitCrit) {
379  if (crit_name == crit->getName()) {
380  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
381  return critKeeper->getValues();
382  }
383  }
384 
385  return em;
386  } //getCriteriaValues
387 
394  std::vector<std::map < std::string , float >> getCriteriaAllValues(std::string crit_name) {
395  std::vector<std::map < std::string , float >> em;
396  if (mSaveCriteriaValues != true) {
397  return em;
398  }
399 
400  for (auto crit : mTwoHitCrit) {
401  if (crit_name == crit->getName()) {
402  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
403  return critKeeper->getAllValues();
404  }
405  }
406 
407  for (auto crit : mThreeHitCrit) {
408  if (crit_name == crit->getName()) {
409  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
410  return critKeeper->getAllValues();
411  }
412  }
413 
414  return em;
415  } // getCriteriaAllValues
416 
423  std::vector<int> getCriteriaTrackIds(std::string crit_name) {
424  std::vector<int> em;
425  if (mSaveCriteriaValues != true) {
426  return em;
427  }
428 
429  for (auto crit : mTwoHitCrit) {
430  if (crit_name == crit->getName()) {
431  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
432  return critKeeper->getTrackIds();
433  }
434  }
435 
436  for (auto crit : mThreeHitCrit) {
437  if (crit_name == crit->getName()) {
438  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
439  return critKeeper->getTrackIds();
440  }
441  }
442 
443  return em;
444  } //getCriteriaTrackIds
445 
451  if (mSaveCriteriaValues != true) {
452  return;
453  }
454 
455  for (auto crit : mTwoHitCrit) {
456  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
457  critKeeper->clear();
458  }
459 
460  for (auto crit : mThreeHitCrit) {
461  auto critKeeper = static_cast<CriteriaKeeper *>(crit);
462  critKeeper->clear();
463  }
464  } // clearSavedCriteria
465 
472  size_t nHitsInHitMap(FwdDataSource::HitMap_t &hitmap) {
473  size_t n = 0;
474  for (auto kv : hitmap) {
475  n += kv.second.size();
476  }
477  return n;
478  }
479 
486  void removeHits(FwdDataSource::HitMap_t &hitmap, std::vector<Seed_t> &tracks) {
487 
488  for (auto track : tracks) {
489  for (auto hit : track) {
490  int sector = hit->getSector();
491 
492  auto hitit = std::find(hitmap[sector].begin(), hitmap[sector].end(), hit);
493 
494  if (hitit != hitmap[sector].end()) {
495  hitmap[sector].erase(hitit);
496  mTotalHitsRemoved++;
497  }
498 
499  } // loop on hits in track
500  } // loop on track
501  } // removeHits
502 
503 
510  void mergeHitmaps( FwdDataSource::HitMap_t &hitmap1, FwdDataSource::HitMap_t &hitmap2 ){
511  static const int numFstSectors = 3;
512  for ( auto kv : hitmap2 ){
513  for ( auto hit : kv.second ){
514  dynamic_cast<FwdHit*>( hit )->setSector( hit->getSector() + numFstSectors );
515  hitmap1[kv.first + numFstSectors].push_back( hit );
516  }
517  }
518  } // mergeHitmaps
519 
523  void cleanup() {
524  /************** Cleanup ****************************************/
525  // Moved cleanup to the start of doEvent, so that the fit results
526  // persist after the call
527  mTrackSeeds.clear();
528  mTrackResults.clear();
529  mEventStats.reset();
530  mTotalHitsRemoved = 0;
531  /************** Cleanup **************************/
532  }
533 
534  bool useMcTrackFinding(){
535  /*************************************************************/
536  // Determine if we should use MC seed finding
537  /*************************************************************/
538  bool mcTrackFinding = true;
539  if (mConfig.exists("TrackFinder")){
540  mcTrackFinding = false;
541  }
542  if (mConfig.exists("TrackFinder") && mConfig.get<bool>( "TrackFinder:mc", false ) == false ){
543  mcTrackFinding = false;
544  }
545  if (mConfig.exists("TrackFinder") && mConfig.get<bool>( "TrackFinder:active", true ) == false){
546  mcTrackFinding = true;
547  }
548  return mcTrackFinding;
549  }
550 
556  void findTrackSeeds() {
557  cleanup();
558  /*************************************************************/
559  // Get the hitmaps
560  /*************************************************************/
561  long long itStart = FwdTrackerUtils::nowNanoSecond();
562  FwdDataSource::HitMap_t &fttHitmap = mDataSource->getFttHits();
563  FwdDataSource::HitMap_t &fstHitmap = mDataSource->getFstHits();
564 
565  /*************************************************************/
566  // Determine seed finding mode
567  /*************************************************************/
568  string hitmapSource = mConfig.get<string>("TrackFinder:source", "ftt");
569  LOG_INFO << "Hitmap Source: " << hitmapSource << endm;
570  mSeedSource = kSeqSeed; // default to FST
571  if (hitmapSource == "fst")
572  mSeedSource = kFstSeed;
573  else if (hitmapSource == "ftt")
574  mSeedSource = kFttSeed;
575  else if (hitmapSource == "seq")
576  mSeedSource = kSeqSeed;
577  else if (hitmapSource == "sim")
578  mSeedSource = kSimSeed;
579  LOG_INFO << "Performing Fwd Seed finding with mode: " << mConfig.get<string>("TrackFinder:source", "ftt") << " = " << mSeedSource << endm;
580  FwdDataSource::McTrackMap_t &mcTrackMap = mDataSource->getMcTracks();
581 
582  long long duration = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
583 
584  mEventStats.mStep1Duration.push_back( duration );
585 
586  /*************************************************************/
587  // DO MC Track Finding (if set to do so)
588  if (useMcTrackFinding()) {
589  doMcTrackFinding(mcTrackMap, mSeedSource);
590  mEventStats.mNumSeeds = mTrackSeeds.size();
591  long long duration2 = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
592  mEventStats.mSeedFindingDuration.push_back( duration2 );
593  return;
594  } else {
595  LOG_DEBUG << "Performing Standard Track Finding" << endm;
596  }
597  /*************************************************************/
598 
599  /*************************************************************/
600  // Standard Track Finding
601  size_t nIterations = mConfig.get<size_t>("TrackFinder:nIterations", 0);
602  for (size_t iIteration = 0; iIteration < nIterations; iIteration++) {
603  if ( mSeedSource == kSimSeed){
604  mergeHitmaps( fstHitmap, fttHitmap );
605  } else {
606  if ( mSeedSource == kFstSeed || mSeedSource == kSeqSeed ){
607  doSeedFindingIteration(iIteration, fstHitmap);
608  }
609  if ( mSeedSource == kFttSeed || mSeedSource == kSeqSeed){
610  doSeedFindingIteration(iIteration, fttHitmap);
611  }
612  }
613  } // iIteration
614  /*************************************************************/
615  mEventStats.mNumSeeds = mTrackSeeds.size();
616  long long duration2 = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
617  mEventStats.mSeedFindingDuration.push_back( duration2 );
618  } // FindTrackSeeds
619 
620 
621  std::vector< genfit::GFRaveVertex * > findFwdVertices( const vector<GenfitTrackResult> &globalTracks ){
622  // we will return a vector of vertices
623  std::vector< genfit::GFRaveVertex * > raveVertices;
624 
625 
626  // The RAVE factory needs the (bare) track pointers for vertex finding
627  vector<genfit::Track*> tracks;
628  for ( auto gtr : globalTracks ){
629  if ( gtr.mTrack ){
630  tracks.push_back( gtr.mTrack.get() );
631  }
632  }
633 
634  bool useBeamConstraint = false;
635  if ( useBeamConstraint ){
636  // TODO: load "official" beamline constraint parameters
637  TMatrixDSym bscm(3);
638  const double bssXY = 2.0;
639  bscm(0, 0) = bssXY*bssXY;
640  bscm(1, 1) = bssXY*bssXY;
641  bscm(2, 2) = 50.5 * 50.5;
642  mGFRVertices.setBeamspot( TVector3( 0, 0, 0 ), bscm );
643  }
644 
645 
646 
647  mGFRVertices.findVertices( &raveVertices, tracks, useBeamConstraint );
648  LOG_DEBUG << "raveVertices.size() = " << raveVertices.size() << endm;
649  for ( auto vert : raveVertices ){
650  LOG_DEBUG << TString::Format( "GFRaveVertex vertex @(%f, %f, %f)\n\n", vert->getPos().X(), vert->getPos().Y(), vert->getPos().Z() ) << endm;
651  LOG_DEBUG << "GFRaveVertex\n";
652  // LOG_DEBUG << "Position: "; vert->getPos().Print();
653  // LOG_DEBUG << "Covariance: "; vert->getCov().Print();
654  LOG_DEBUG << "Ndf: " << vert->getNdf() << ", Chi2: " << vert->getChi2() << ", Id: " << vert->getId() << "\n";
655  LOG_DEBUG << "Number of tracks: " << vert->getNTracks() << "\n";
656  }
657  // if there is no Event Vertex, then we will use the first vertex found
658  // if ( raveVertices.size() > 0 ){
659  // mEventVertex = raveVertices[0]->getPos();
660  // // mEventVertexHit.setPos( mEventVertex );
661  // }
662 
663 
664 
665  return raveVertices;
666  }
667 
675  GenfitTrackResult fitTrack(Seed_t &seed, TVector3 *momentumSeedState = nullptr) {
676  LOG_DEBUG << "FwdTracker::fitTrack->" << endm;
677  mEventStats.mAttemptedFits++;
678  // We will build this up as we go
679  GenfitTrackResult gtr;
680 
681  LOG_DEBUG << "--Setting seed on GenfitTrackResult, seed has " << seed.size() << " hits" << endm;
682  // First, set the seed information
683  gtr.setSeed( seed );
684 
685  // If we are using a provided momentum state
686  if ( momentumSeedState ){
687  LOG_DEBUG << "--FitTrack with provided momentum seed state" << endm;
688  mTrackFitter->fitTrack( seed, momentumSeedState );
689  } else {
690  LOG_DEBUG << "--FitTrack without provided momentum seed state" << endm;
691  mTrackFitter->fitTrack( seed );
692  }
693 
694  /*******************************************************/
695  // Get the track from the fitter
696  // and set the track in the GenfitTrackResult
697  if (mTrackFitter->getTrack() != nullptr ){
698  LOG_DEBUG << "--FitTrack found, setting seed and track only" << endm;
699  gtr.set( seed, mTrackFitter->getTrack() );
700 
701  if (gtr.mStatus && gtr.mStatus->isFitConvergedFully()) {
702  mEventStats.mGoodFits++;
703  } else {
704  mEventStats.mFailedFits++;
705  }
706  } else { // set the track as a failed fit, but keep the seed info
707  LOG_DEBUG << "--FitTrack failed, setting seed only" << endm;
708  mEventStats.mFailedFits++;
709  }
710  LOG_DEBUG << "<-FwdTracker::fitTrack complete" << endm;
711  return gtr;
712  } // fitTrack
713 
714 
715 
728  void doTrackFitting( const std::vector<Seed_t> &trackSeeds) {
729  LOG_DEBUG << ">>doTrackFitting" << endm;
730  if (!mDoTrackFitting)
731  return;
732  long long itStart = FwdTrackerUtils::nowNanoSecond();
733 
734  std::vector<GenfitTrackResult> globalTracks;
735  std::vector<GenfitTrackResult> primaryTracks;
736 
737  // Should we try to refit the track with aadditional points from other detectors?
738  const bool doRefit = mConfig.get<bool>("TrackFitter:refit", false);
739  LOG_INFO << "TrackFitter:refit = " << doRefit << endm;
740 
741  // Should we use the MC momentum as a seed for the fit?
742  const bool useMcSeedMomentum = mConfig.get<bool>("TrackFitter:mcSeed", false);
743  LOG_INFO << "TrackFitter:mcSeed = " << useMcSeedMomentum << endm;
744 
745  LOG_DEBUG << "Starting track fitting loop, mTrackResults.size() = " << mTrackResults.size() << endm;
746  LOG_DEBUG << "Starting Track fitting loop on " << trackSeeds.size() << " track seeds" << endm;
747  size_t index = 0;
748  for (auto t : trackSeeds) {
749 
750  GenfitTrackResult gtrGlobalRefit; // will store refit if needed
751  LOG_DEBUG << "\tTrack seed initial global fit #" << index << endm;
752  /***********************************************************************************************************/
753  // Tracking Step 1
754  // Fit each accepted track seed
755 
756  // If we are using MC momentum get it from associated track
757  TVector3 momentumSeedStateMc;
758  TVector3 *pMomSeedState = nullptr;
759  int idt = 0;
760  double qual = 0;
761  // Get the quality and MC truth id
762  idt = MCTruthUtils::dominantContribution(t, qual);
763  LOG_INFO << "\t\tMc Match idTruth=" << idt << ", quality = " << qual << endm;
764  if (true == useMcSeedMomentum) {
765  /*******************************************************/
766  // Only for Simulation
767  // Calculate the MC info first and check filters
768 
769  auto mctm = mDataSource->getMcTracks();
770  // get the MC track momentum if we can (may be used for state seed)
771  if (mctm.count(idt)) {
772  auto mct = mctm[idt];
773  momentumSeedStateMc.SetPtEtaPhi(mct->mPt, mct->mEta, mct->mPhi);
774  // pMomSeedState = &momentumSeedStateMc;
775  } else {
776  LOG_WARN << "\tRequested MC momentum for fit seed, but couldnt find MC Track for id: " << idt << ", qual: " << qual << endm;
777  }
778  /*******************************************************/
779  } // else if pMomSeedState = nullptr, a seed momentum will be computed from seed points by tracker
780 
781  // Fit the track seed and get the GenfitTrackResult
782  GenfitTrackResult gtrGlobal = fitTrack(t, pMomSeedState);
783  gtrGlobal.setDCA( mEventVertex );
784 
785  LOG_DEBUG << "\tFit track seed with " << gtrGlobal.mSeed.size() << " hits" << endm;
786  LOG_DEBUG << "\t\t McTrack Id = " << gtrGlobal.mIdTruth << ", QA = " << gtrGlobal.mQaTruth << endm;
787  // End Step 1
788  /*******************************************************/
789 
790 
791  // if the first fit fails then we cannot proceed with the refit steps
792  if (gtrGlobal.mIsFitConvergedPartially == false) {
793  LOG_WARN << "\tInitial fitting failed for seed " << index << endm;
794  LOG_DEBUG << "\tFitting failed for seed " << index << endm;
795  LOG_DEBUG << "\tSkipping the refit steps but saving the seed and failed fit" << endm;
796  globalTracks.push_back( gtrGlobal );
797  index++;
798  continue;
799  // BREAK OUT OF THE LOOP
800  }
801  if (doRefit == false) {
802  LOG_INFO << "\tRefit is disabled, saving the seed and initial fit" << endm;
803  gtrGlobal.isPrimary = false;
804  globalTracks.push_back( gtrGlobal );
805  index++;
806  continue;
807  // BREAK OUT OF THE LOOP
808  }
809 
810  /***********************************************************************************************************/
811  // Tracking Step 2
812  // Look for additional hits in the other tracking detector
813  // and add the new hits to the track
814 
815  // If requested, use the MC track finding to add hits to the track
816  if (useMcTrackFinding()) {
817  if (mSeedSource != kFttSeed) addFttHitsMc( gtrGlobal );
818  if (mSeedSource != kFstSeed) addFstHitsMc( gtrGlobal );
819  // globalTracks.push_back( gtrGlobalRefit );
820  // index++;
821  // continue; // below is for "real" track finding only, for MC jump to next track
822  } else {
823  // If we are not using MC track finding,
824  // then we will look for additional hits via projections
825  if (mSeedSource != kFttSeed){ // Look for FTT hits if it was not the original seed source
826  for ( int i = 0; i < FwdSystem::sNFttLayers; i++ ){
827  addFttHits( gtrGlobal, i );
828  }
829  }
830  if (mSeedSource != kFstSeed ){ // Look for FST hits if it was not the original seed source
831  for ( int i = 0; i < FwdSystem::sNFstLayers; i++ ){
832  addFstHits( gtrGlobal, i );
833  }
834  }
835  // global refit
836  }
837  // End Step 2
838  /***********************************************************************************************************/
839 
840  /***********************************************************************************************************/
841  // Tracking Step 3: Fit the new global track with additional hits
842  gtrGlobalRefit = fitTrack( gtrGlobal.mSeed, &gtrGlobal.mMomentum );
843  gtrGlobalRefit.setDCA( mEventVertex );
844  // End Step 3
845  /***********************************************************************************************************/
846 
847  /***********************************************************************************************************/
848  // Tracking Step 4: Save the best global track result
849  GenfitTrackResult *activeTrack = &gtrGlobal;
850  if ( gtrGlobalRefit.mIsFitConvergedPartially ){
851  activeTrack = &gtrGlobalRefit;
852  }
853 
854  if ( !activeTrack->mIsFitConvergedPartially ){
855  // should not be possible according to above logic...
856  LOG_WARN << "\tFWD global track fit failed (both initial + refit)" << endm;
857  continue;
858  }
859 
860  // if the refit is successful,
861  // then we will add the track to the globalTracks
862  // if not keep the original global track
863  activeTrack->isPrimary = false; // should be default but just make sure
864  globalTracks.push_back( *activeTrack ); // save this global track result
865  // End Step 4
866  /***********************************************************************************************************/
867  } // loop on track seeds
868 
869  long long duration1 = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
870  mEventStats.mFitDuration.push_back( duration1 );
871  float perGood = (float) mEventStats.mGoodFits / (float) mEventStats.mAttemptedFits;
872  float perFailed = (float) mEventStats.mFailedFits / (float) mEventStats.mAttemptedFits;
873  float perCardinals = (float) mEventStats.mGoodCardinals / (float) mEventStats.mGoodFits;
874  LOG_DEBUG << "\tGlobal Track Fitting Results:" <<
875  TString::Format(
876  "Attempts = %d, Good = %d (%f%%), Failed = %d (%f%%), GoodCardinals = %d (%f%%)",
877  mEventStats.mAttemptedFits,
878  mEventStats.mGoodFits,
879  perGood,
880  mEventStats.mFailedFits,
881  perFailed,
882  mEventStats.mGoodCardinals,
883  perCardinals ) << endm;
884 
885  const bool do_fwd_vertex_finding = true;
886  if (do_fwd_vertex_finding){
887  /***********************************************************************************************************/
888  // Step 5: Find the FWD Vertices
889  LOG_DEBUG << "\tStarting Track Fitting Step 3 (FWD Vertex Finding)" << endm;
890  auto fwdVertices = findFwdVertices( globalTracks );
891  mEventStats.mNumFwdVertices = fwdVertices.size();
892 
893  for ( auto vert : fwdVertices ){
894  LOG_DEBUG << "\tFound FWD Vertex @(" << vert->getPos().X() << ", " << vert->getPos().Y() << ", " << vert->getPos().Z() << ")" << endm;
895  LOG_DEBUG << "\t\tvs mEventVertexHit: " << mEventVertexHit.getX() << ", " << mEventVertexHit.getY() << ", " << mEventVertexHit.getZ() << endm;
896  }
897 
898  // End Step 5
899  /***********************************************************************************************************/
900  } else {
901  LOG_INFO << "Event configuration is skipping FWD vertex finding" << endm;
902  }
903 
904 
905  const bool do_fwd_primary_fitting = true;
906  /***********************************************************************************************************/
907  // Step 6: Refit the track with the primary vertex
908  index = 0;
909  if (do_fwd_primary_fitting){
910  // Now try refitting every track with the primary vertex
911  for (auto &gtr : globalTracks) {
912  LOG_INFO << "Refitting Track " << index << ", McId=" << gtr.mIdTruth << " with Primary Vertex, seed already has: " << gtr.mSeed.size() << " hits" << endm;
913  LOG_INFO << "mEventVertexHit: " << mEventVertexHit.getX() << ", " << mEventVertexHit.getY() << ", " << mEventVertexHit.getZ() << endm;
914  // just use the global track to build the track that will use the PV also
915  Seed_t seedWithPV = gtr.mSeed;
916  seedWithPV.push_back( &mEventVertexHit );
917 
918  // If we are using MC momentum get it from associated track
919  TVector3 momentumSeedStateMc;
920  TVector3 *pMomSeedState = &gtr.mMomentum;
921  // if (true == useMcSeedMomentum) {
922  // /*******************************************************/
923  // // Only for Simulation
924  // // get the MC track momentum if we can (may be used for state seed)
925  // auto mctm = mDataSource->getMcTracks();
926  // if (mctm.count(gtr.mIdTruth)) {
927  // auto mct = mctm[gtr.mIdTruth];
928  // momentumSeedStateMc.SetPtEtaPhi(mct->mPt, mct->mEta, mct->mPhi);
929  // pMomSeedState = &momentumSeedStateMc;
930  // LOG_INFO << "Setting momentum to MC Seed state value for primary refit" << endm;
931  // } else {}
932  // /*******************************************************/
933  // } // else if pMomSeedState = nullptr, a seed momentum will be computed from seed points by tracker
934 
935  GenfitTrackResult gtrPV = fitTrack(seedWithPV, pMomSeedState);
936  if ( gtrPV.mIsFitConvergedFully ){
937  mEventStats.mGoodPrimaryFits++;
938  } else {
939  mEventStats.mFailedPrimaryFits++;
940  continue;
941  }
942  gtrPV.setDCA( mEventVertex );
943  gtrPV.isPrimary = true;
944  mEventStats.mAttemptedPrimaryFits ++;
945  primaryTracks.push_back( gtrPV );
946  index++;
947  }
948  // End Step 6
949  /***********************************************************************************************************/
950  } else {
951  LOG_INFO << "Event configuration is skipping primary track fitting" << endm;
952  }
953 
954  float perGoodPrim = (float) mEventStats.mGoodPrimaryFits / (float) mEventStats.mAttemptedPrimaryFits;
955  float perFailedPrim = (float) mEventStats.mFailedPrimaryFits / (float) mEventStats.mAttemptedPrimaryFits;
956  LOG_DEBUG << "\tPrimary Track Fitting Results:" <<
957  TString::Format(
958  "Attempts = %d, Good = %d (%f%%), Failed = %d (%f%%)",
959  mEventStats.mAttemptedPrimaryFits,
960  mEventStats.mGoodPrimaryFits,
961  perGoodPrim,
962  mEventStats.mFailedPrimaryFits,
963  perFailedPrim
964  ) << endm;
965 
966  // Add the global and primary tracks to the results
967  LOG_DEBUG << "Ending track fitting loop, mTrackResults.size() = " << mTrackResults.size() << endm;
968  mTrackResults.insert( mTrackResults.end(), globalTracks.begin(), globalTracks.end() );
969  LOG_DEBUG << "Copied globals, now mTrackResults.size() = " << mTrackResults.size() << endm;
970  mTrackResults.insert( mTrackResults.end(), primaryTracks.begin(), primaryTracks.end() );
971 
972  long long duration2 = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
973  LOG_DEBUG << "Track fitting took " << duration2 << "ms" << endm;
974  LOG_DEBUG << "We fit " << globalTracks.size() << " global tracks and " << primaryTracks.size() << " primary tracks, total = " << mTrackResults.size() << endm;
975  } // doTrackFitting
976 
983  void doMcTrackFinding(FwdDataSource::McTrackMap_t &mcTrackMap, int seedSource) {
984  LOG_INFO << "Running MC Seed Finding, mode: " << seedSource << endm;
985 
986  // If we want sequential MC track finding then do them each individually
987  if ( seedSource == kSeqSeed ){
988  doMcTrackFinding( mcTrackMap, kFstSeed );
989  doMcTrackFinding( mcTrackMap, kFttSeed );
990  return;
991  }
992 
993  mTrackSeedsThisIteration.clear();
994  // we will build reco tracks from each McTrack
995  for (auto kv : mcTrackMap) {
996 
997  auto mc_track = kv.second;
998  LOG_DEBUG << "McTrack[ " << kv.first << " ]: nFtt=" << mc_track->mFttHits.size() << ", nFst=" << mc_track->mFstHits.size() << endm;
999 
1000  if (seedSource == kFttSeed && mc_track->mFttHits.size() < 2){ // require min 4 FTT hits on track
1001  continue;
1002  }
1003 
1004  if (seedSource == kFstSeed && mc_track->mFstHits.size() < 2 ) { // require min 3 FST hits on track
1005  LOG_DEBUG << "Skipping McSeedFinding bc FST hits < 2" << endm;
1006  continue;
1007  }
1008 
1009  if ( seedSource == kSimSeed && mc_track->mFstHits.size() < 2 && mc_track->mFttHits.size() < 2 ){
1010  continue;
1011  }
1012 
1013  std::set<size_t> uvid;
1014  Seed_t track;
1015 
1016  if ( seedSource != kFstSeed ){ // FTT is used unless we are ONLY considering FST
1017  for (auto h : mc_track->mFttHits) {
1018  track.push_back(h);
1019  uvid.insert(static_cast<FwdHit *>(h)->_vid);
1020  }
1021  }
1022  if (seedSource != kFttSeed ) { // FST
1023  for (auto h : mc_track->mFstHits) {
1024  track.push_back(h);
1025  uvid.insert(static_cast<FwdHit *>(h)->_vid);
1026  }
1027  }
1028 
1029  if (uvid.size() == track.size()) { // only add tracks that have one hit per volume
1030  mTrackSeedsThisIteration.push_back(track);
1031  int idt = 0;
1032  double qual = 0;
1033  idt = MCTruthUtils::dominantContribution(track, qual);
1034  } else {
1035  //Skipping track that doesnt have hits on all layers
1036  }
1037  }
1038 
1039  LOG_DEBUG << "McTrackFinding Found: " << mTrackSeedsThisIteration.size() << " tracks" << endm;
1040  // doTrackFitting(mTrackSeedsThisIteration);
1041 
1042  // Now save to the main reco track list
1043  mTrackSeeds.insert( mTrackSeeds.end(), mTrackSeedsThisIteration.begin(), mTrackSeedsThisIteration.end() );
1044  } //doMcTrackFinding
1045 
1056  size_t sliceHitMapInPhi( FwdDataSource::HitMap_t &inputMap, FwdDataSource::HitMap_t &outputMap, float phi_min, float phi_max ){
1057  size_t n_hits_kept = 0;
1058 
1059  outputMap.clear(); // child STL containers will get cleared too
1060  for ( auto kv : inputMap ){
1061  for ( KiTrack::IHit* hit : kv.second ){
1062  TVector3 vec(hit->getX(), hit->getY(), hit->getZ() );
1063  if ( vec.Phi() < phi_min || vec.Phi() > phi_max ) continue;
1064 
1065  // now add the hits to the sliced map
1066  outputMap[kv.first].push_back( hit );
1067  n_hits_kept ++;
1068  } // loop on hits
1069  } // loop on map
1070  return n_hits_kept;
1071  } // sliceHitMapInPhi
1072 
1079  vector<Seed_t> doSeedFindingOnHitmapSubset( size_t iIteration, FwdDataSource::HitMap_t &hitmap ) {
1080  long long itStart = FwdTrackerUtils::nowNanoSecond();
1081 
1082  std::vector<Seed_t> acceptedTrackSeeds;
1083  /*************************************************************/
1084  // Step 2
1085  // build 2-hit segments (setup parent child relationships)
1086  /*************************************************************/
1087  // Initialize the segment builder with sorted hits
1088  KiTrack::SegmentBuilder builder(hitmap);
1089 
1090  // Load the criteria used for 2-hit segments
1091  // This loads from XML config if available
1092  std::string criteriaPath = "TrackFinder.Iteration[" + std::to_string(iIteration) + "].SegmentBuilder";
1093 
1094  if (false == mConfig.exists(criteriaPath)) {
1095  // Use the default for all iterations if it is given.
1096  // If not then no criteria will be applied
1097  criteriaPath = "TrackFinder.SegmentBuilder";
1098  }
1099 
1100  clearCriteria( mTwoHitCrit );
1101  mTwoHitCrit = loadCriteria(criteriaPath);
1102  builder.addCriteria(mTwoHitCrit);
1103 
1104  // Setup the connector (this tells it how to connect hits together into segments)
1105  std::string connPath = "TrackFinder.Iteration[" + std::to_string(iIteration) + "].Connector";
1106 
1107  if (false == mConfig.exists(connPath))
1108  connPath = "TrackFinder.Connector";
1109 
1110  unsigned int distance = mConfig.get<unsigned int>(connPath + ":distance", 1);
1111  if (mSeedSource == kFttSeed){
1112  distance = 2; // set distance to 2 for FTT
1113  }
1114 
1115  FwdConnector connector(distance);
1116  builder.addSectorConnector(&connector);
1117  LOG_DEBUG << "Connector added: " << endm;
1118  // Get the segments and return an automaton object for further work
1119 
1120  KiTrack::Automaton automaton = builder.get1SegAutomaton();
1121  LOG_DEBUG << TString::Format( "nSegments=%lu", automaton.getSegments().size() ).Data() << endm;
1122  LOG_DEBUG << TString::Format( "nConnections=%u", automaton.getNumberOfConnections() ).Data() << endm;
1123 
1124  if (automaton.getNumberOfConnections() > 9000 ){
1125  LOG_ERROR << "Got too many connections, bailing out of tracking" << endm;
1126  return acceptedTrackSeeds;
1127  }
1128 
1129  // at any point we can get a list of tracks out like this:
1130  // std::vector < std::vector< KiTrack::IHit* > > tracks = automaton.getTracks();
1131  // we can apply an optional parameter <nHits> to only get tracks with >=nHits in them
1132 
1133  long long duration = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
1134  mEventStats.mStep2Duration.push_back(duration);
1135  itStart = FwdTrackerUtils::nowNanoSecond();
1136 
1137  /*************************************************************/
1138  // Step 3
1139  // build 3-hit segments from the 2-hit segments
1140  /*************************************************************/
1141  automaton.clearCriteria();
1142  automaton.resetStates();
1143  criteriaPath = "TrackFinder.Iteration[" + std::to_string(iIteration) + "].ThreeHitSegments";
1144 
1145  if (false == mConfig.exists(criteriaPath))
1146  criteriaPath = "TrackFinder.ThreeHitSegments";
1147 
1148  clearCriteria( mThreeHitCrit );
1149  mThreeHitCrit = loadCriteria(criteriaPath);
1150  automaton.addCriteria(mThreeHitCrit);
1151 
1152  duration = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
1153  mEventStats.mStep3Duration.push_back( duration );
1154  if (duration > 2000 || automaton.getNumberOfConnections() > 9000){
1155  LOG_WARN << "The Three Hit Criteria took more than 200ms to process, duration: " << duration << " ms" << endm;
1156  LOG_WARN << "bailing out (skipping subset HNN)" << endm;
1157 
1158  std::string subsetPath = "TrackFinder.Iteration[" + std::to_string(iIteration) + "].SubsetNN";
1159  size_t minHitsOnTrack = mConfig.get<size_t>(subsetPath + ":min-hits-on-track", FwdSystem::sNFttLayers);
1160  acceptedTrackSeeds = automaton.getTracks(minHitsOnTrack);
1161  return acceptedTrackSeeds;
1162  }
1163  itStart = FwdTrackerUtils::nowNanoSecond();
1164 
1165  LOG_DEBUG << TString::Format( "nSegments=%lu", automaton.getSegments().size() ).Data() << endm;
1166  LOG_DEBUG << TString::Format( "nConnections=%u", automaton.getNumberOfConnections() ).Data() << endm;
1167 
1168  /*************************************************************/
1169  // Step 4
1170  // Get the tracks from the possible tracks that are the best subset
1171  /*************************************************************/
1172  std::string subsetPath = "TrackFinder.Iteration[" + std::to_string(iIteration) + "].SubsetNN";
1173 
1174  if (false == mConfig.exists(subsetPath))
1175  subsetPath = "TrackFinder.SubsetNN";
1176 
1177  // only for debug really
1178  bool findSubsets = mConfig.get<bool>(subsetPath + ":active", true);
1179 
1180  if (findSubsets) {
1181  size_t minHitsOnTrack = mConfig.get<size_t>(subsetPath + ":min-hits-on-track", 7);
1182  // Getting all tracks with at least minHitsOnTrack hits on them
1183  std::vector<Seed_t> tracks = automaton.getTracks(minHitsOnTrack);
1184 
1185  float omega = mConfig.get<float>(subsetPath + ".Omega", 0.75);
1186  float stableThreshold = mConfig.get<float>(subsetPath + ".StableThreshold", 0.1);
1187  float Ti = mConfig.get<float>(subsetPath + ".InitialTemp", 2.1);
1188  float Tf = mConfig.get<float>(subsetPath + ".InfTemp", 0.1);
1189 
1190  KiTrack::SubsetHopfieldNN<Seed_t> subset;
1191  subset.add(tracks);
1192  subset.setOmega(omega);
1193  subset.setLimitForStable(stableThreshold);
1194  subset.setTStart(Ti);
1195  subset.setTInf(Tf);
1196 
1197  SeedCompare comparer;
1198  SeedQual quality;
1199 
1200  subset.calculateBestSet(comparer, quality);
1201 
1202  acceptedTrackSeeds = subset.getAccepted();
1203 
1204  // this call takes a long time due to possible huge combinatorics.
1205  // rejectedTracks = subset.getRejected();
1206  // LOG_DEBUG << "We had " << tracks.size() << " tracks. Accepted = " << acceptedTrackSeeds.size() << ", Rejected = " << rejectedTracks.size() << endm;
1207 
1208  } else { // the subset and hit removal
1209  size_t minHitsOnTrack = mConfig.get<size_t>(subsetPath + ":min-hits-on-track", FwdSystem::sNFttLayers);
1210  acceptedTrackSeeds = automaton.getTracks(minHitsOnTrack);
1211  }// subset off
1212 
1213  duration = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
1214  mEventStats.mStep4Duration.push_back( duration );
1215  if (duration > 500){
1216  LOG_WARN << "The took more than 500ms to process, duration: " << duration << " ms" << endm;
1217  LOG_WARN << "We got " << acceptedTrackSeeds.size() << " tracks this round" << endm;
1218  }
1219  LOG_DEBUG << "We got " << acceptedTrackSeeds.size() << " tracks this round" << endm;
1220  return acceptedTrackSeeds;
1221  } // doSeedFindingOnHitmapSubset
1222 
1229  void doSeedFindingIteration(size_t iIteration, FwdDataSource::HitMap_t &hitmap) {
1230 
1231  // empty the list of reco tracks for the iteration
1232  mTrackSeedsThisIteration.clear();
1233 
1234  // check to see if we have hits!
1235  size_t nHitsThisIteration = nHitsInHitMap(hitmap);
1236 
1237  const int minHitsToConsider = 3;
1238  if (nHitsThisIteration < minHitsToConsider) {
1239  // No hits left in the hitmap! Skipping this iteration
1240  LOG_INFO << "No hits to consider in this iteration, skipping" << endm;
1241  return;
1242  }
1243 
1244  std::string pslPath = "TrackFinder.Iteration["+ std::to_string(iIteration) + "]:nPhiSlices";
1245  if ( false == mConfig.exists( pslPath ) ) pslPath = "TrackFinder:nPhiSlices";
1246  size_t phi_slice_count = mConfig.get<size_t>( pslPath, 1 );
1247 
1248  if ( phi_slice_count <= 1 || phi_slice_count > 100 ) { // no phi slicing!
1249  LOG_DEBUG << "Tracking without phi slices" << endm;
1250  /*************************************************************/
1251  // Steps 2 - 4 here
1252  /*************************************************************/
1253  auto acceptedTracks = doSeedFindingOnHitmapSubset( iIteration, hitmap );
1254  mTrackSeedsThisIteration.insert( mTrackSeedsThisIteration.end(), acceptedTracks.begin(), acceptedTracks.end() );
1255  } else {
1256 
1257  FwdDataSource::HitMap_t slicedHitMap;
1258 
1259  if ( phi_slice_count == 0 || phi_slice_count > 100 ){
1260  LOG_WARN << "Invalid phi_slice_count = " << phi_slice_count << ", resetting to 1" << endm;
1261  phi_slice_count= 1;
1262  }
1263  float phi_slice = 2 * TMath::Pi() / (float) phi_slice_count;
1264  for ( size_t phi_slice_index = 0; phi_slice_index < phi_slice_count; phi_slice_index++ ){
1265 
1266  float phi_min = phi_slice_index * phi_slice - TMath::Pi();
1267  float phi_max = (phi_slice_index + 1) * phi_slice - TMath::Pi();
1268  LOG_INFO << TString::Format( "phi slice = (%f, %f)", phi_min, phi_max ) << endm;
1269 
1270  /*************************************************************/
1271  // Step 1
1272  // Slice the hitmap into a phi section if needed
1273  // If we do that, check again that we arent wasting time on empty sections
1274  /*************************************************************/
1275  size_t nHitsThisSlice = 0;
1276  if ( phi_slice_count > 1 ){
1277  nHitsThisSlice = sliceHitMapInPhi( hitmap, slicedHitMap, phi_min, phi_max );
1278  if ( nHitsThisSlice < minHitsToConsider ) {
1279  continue;
1280  }
1281  } else { // no need to slice
1282  // I think this incurs a copy, maybe we can find a way to avoid.
1283  slicedHitMap = hitmap;
1284  }
1285 
1286  /*************************************************************/
1287  // Steps 2 - 4 here
1288  /*************************************************************/
1289  auto acceptedTracks = doSeedFindingOnHitmapSubset( iIteration, slicedHitMap );
1290  mTrackSeedsThisIteration.insert( mTrackSeedsThisIteration.end(), acceptedTracks.begin(), acceptedTracks.end() );
1291  } //loop on phi slices
1292  }// if loop on phi slices
1293 
1294  /*************************************************************/
1295  // Step 5
1296  // Remove the hits from any track that was found
1297  /*************************************************************/
1298  std::string hrmPath = "TrackFinder.Iteration["+ std::to_string(iIteration) + "].HitRemover";
1299  if ( false == mConfig.exists( hrmPath ) ) hrmPath = "TrackFinder.HitRemover";
1300 
1301  if ( true == mConfig.get<bool>( hrmPath + ":active", true ) ){
1302  removeHits( hitmap, mTrackSeedsThisIteration );
1303  }
1304 
1305  LOG_DEBUG << " Found " << mTrackSeedsThisIteration.size() << " seed tracks this iteration" << endm;
1306  // Add the set of all accepted tracks (this iteration) to our collection of found tracks from all iterations
1307  mTrackSeeds.insert( mTrackSeeds.end(), mTrackSeedsThisIteration.begin(), mTrackSeedsThisIteration.end() );
1308  } // doSeedFindingIteration
1309 
1315  FwdDataSource::HitMap_t hitmap = mDataSource->getFstHits();
1316  if ( gtr.mStatus->isFitConverged() == false || gtr.mMomentum.Perp() < 1e-3) {
1317  LOG_DEBUG << "Skipping addFstHitsMc, fit failed" << endm;
1318  return;
1319  }
1320  mEventStats.mPossibleReFit++;
1321  Seed_t fstHitsThisTrack;
1322 
1323  for (size_t j = 0; j < 3; j++) {
1324  for (auto h0 : hitmap[j]) {
1325  if (dynamic_cast<FwdHit *>(h0)->_tid == gtr.mIdTruth) {
1326  fstHitsThisTrack.push_back(h0);
1327  break;
1328  }
1329  } // loop on hits in this layer of hitmap
1330  } // loop on hitmap layers
1331 
1332  LOG_DEBUG << "Found " << gtr.mSeed.size() << " existing seed points" << endm;
1333  LOG_DEBUG << "Adding " << fstHitsThisTrack.size() << " new FST seed points" << endm;
1334 
1335  if (fstHitsThisTrack.size() >= 1) {
1336  mEventStats.mAttemptedReFits++;
1337  gtr.mSeed.insert( gtr.mSeed.end(), fstHitsThisTrack.begin(), fstHitsThisTrack.end() );
1338  } // we have 3 Si hits to refit with
1339  } // addFstHitsMc
1340 
1348  void addFttHits( GenfitTrackResult &gtr, size_t disk ) {
1349  FwdDataSource::HitMap_t hitmap = mDataSource->getFttHits();
1350  if ( disk > 3 ) {
1351  LOG_WARN << "Invalid FTT disk number: " << disk << ", cannot add Ftt points to track" << endm;
1352  return;
1353  }
1354  if (gtr.mIsFitConverged != true)
1355  return;
1356  mEventStats.mPossibleReFit++;
1357 
1358  Seed_t hits_near_plane;
1359  try {
1360  auto msp = mTrackFitter->projectToFtt(disk, gtr.mTrack);
1361 
1362  // now look for Ftt hits near the specified state
1363  hits_near_plane = findFttHitsNearProjectedState(hitmap[disk], msp);
1364  LOG_DEBUG << " Found #FTT hits on plane #" << disk << TString::Format( " = [%ld]", hits_near_plane.size() ) << endm;
1365  } catch (genfit::Exception &e) {
1366  // Failed to project
1367  LOG_WARN << "Unable to get Ftt projections: " << e.what() << endm;
1368  }
1369 
1370  LOG_DEBUG << "Found " << gtr.mSeed.size() << " existing seed points" << endm;
1371 
1372  if ( hits_near_plane.size() > 0 ){
1373  mEventStats.mAttemptedReFits++;
1374  LOG_DEBUG << "Adding " << hits_near_plane.size() << " new FTT seed points" << endm;
1375  gtr.mSeed.insert( gtr.mSeed.end(), hits_near_plane.begin(), hits_near_plane.end() );
1376  }
1377  return;
1378  } // addFttHits
1379 
1385  LOG_DEBUG << "Looking for FTT hits on this track (MC lookup)" << endm;
1386  LOG_DEBUG << "Track TruthId = " << gtr.mIdTruth << " vs. " << gtr.mTrack->getMcTrackId() << endm;
1387  FwdDataSource::HitMap_t hitmap = mDataSource->getFttHits();
1388  if ( gtr.mStatus->isFitConverged() == false || gtr.mMomentum.Perp() < 1e-6) {
1389  LOG_DEBUG << "Skipping addFttHitsMc on this track, fit failed" << endm;
1390  return;
1391  }
1392 
1393  mEventStats.mPossibleReFit++;
1394  Seed_t fttHitsForThisTrack;
1395 
1396  for (size_t j = 0; j < 4; j++) {
1397  for (auto h0 : hitmap[j]) {
1398  if (dynamic_cast<FwdHit *>(h0)->_tid == gtr.mIdTruth) {
1399  fttHitsForThisTrack.push_back( h0 );
1400  break;
1401  }
1402  } // loop on hits in this layer of hitmap
1403  } // loop on hitmap layers
1404 
1405  LOG_DEBUG << "Found " << fttHitsForThisTrack.size() << " FTT Hits on this track (MC lookup)" << endm;
1406 
1407  if (fttHitsForThisTrack.size() >= 1) {
1408  mEventStats.mAttemptedReFits++;
1409  gtr.mSeed.insert( gtr.mSeed.end(), fttHitsForThisTrack.begin(), fttHitsForThisTrack.end() );
1410  } // we have at least one Fst hit to refit with
1411  } // add Ftt hits via MC associations
1412 
1419  void addFstHits( GenfitTrackResult &gtr, size_t disk ) {
1420  FwdDataSource::HitMap_t hitmap = mDataSource->getFstHits();
1421  if (gtr.mIsFitConverged == false) {
1422  // Original Track fit did not converge, skipping
1423  return;
1424  }
1425  if ( disk > 2 ){
1426  LOG_ERROR << "Invalid FST disk number: " << disk << endm;
1427  return;
1428  }
1429  mEventStats.mPossibleReFit++;
1430 
1431  Seed_t nearby_hits;
1432  try {
1433  // get measured state on plane at specified disk
1434  auto msp = mTrackFitter->projectToFst(disk, gtr.mTrack);
1435  // now look for Si hits near this state
1436  nearby_hits = findFstHitsNearProjectedState(hitmap[disk], msp);
1437  } catch (genfit::Exception &e) {
1438  LOG_WARN << "Unable to get projections: " << e.what() << endm;
1439  }
1440  LOG_DEBUG << "Track already has " << gtr.mSeed.size() << " existing seed points" << endm;
1441 
1442  if ( nearby_hits.size() > 0 ){
1443  mEventStats.mAttemptedReFits++;
1444  LOG_DEBUG << "Adding " << nearby_hits.size() << " new FST seed points from disk " << disk << endm;
1445  gtr.mSeed.insert( gtr.mSeed.end(), nearby_hits.begin(), nearby_hits.end() );
1446  }
1447  return;
1448  } // addFstHits
1449 
1459  Seed_t findFstHitsNearProjectedState(Seed_t &available_hits, genfit::MeasuredStateOnPlane &msp, double dphi = 0.004 * 20.5, double dr = 2.75 * 2) {
1460  double probe_phi = TMath::ATan2(msp.getPos().Y(), msp.getPos().X());
1461  double probe_r = sqrt(pow(msp.getPos().X(), 2) + pow(msp.getPos().Y(), 2));
1462 
1463  Seed_t found_hits;
1464 
1465  for (auto h : available_hits) {
1466  double h_phi = TMath::ATan2(h->getY(), h->getX());
1467  double h_r = sqrt(pow(h->getX(), 2) + pow(h->getY(), 2));
1468  double mdphi = fabs(h_phi - probe_phi);
1469  if (mdphi > 2*3.1415926)
1470  mdphi = mdphi - 2*3.1415926;
1471 
1472  if ( mdphi < dphi && fabs( h_r - probe_r ) < dr) { // handle 2pi edge
1473  found_hits.push_back(h);
1474  }
1475  }
1476 
1477  return found_hits;
1478  } // findFstHitsNearProjectedState
1479 
1490  Seed_t findFttHitsNearProjectedState(Seed_t &available_hits, genfit::MeasuredStateOnPlane &msp, double dx = 1.5, double dy = 1.5) {
1491 
1492  Seed_t found_hits;
1493  TLorentzVector lv1, lv2;
1494  lv1.SetPxPyPzE( msp.getPos().X(), msp.getPos().Y(), 0, 1 );
1495 
1496  double mindx = 99;
1497  double mindy = 99;
1498  double mindr = 99;
1499  double mindp = 99;
1500  KiTrack::IHit *closest = nullptr;
1501 
1502  for (auto h : available_hits) {
1503 
1504  lv2.SetPxPyPzE( h->getX(), h->getY(), 0, 1 );
1505  double sr = lv1.Pt() - lv2.Pt();
1506  double sp = lv1.DeltaPhi( lv2 );
1507  double sx = h->getX() - msp.getPos().X();
1508  double sy = h->getY() - msp.getPos().Y();
1509 
1510  if ( fabs(sr) < fabs(mindr) )
1511  mindr = sr;
1512  if ( fabs(sp) < fabs(mindp) ){
1513  mindp = sp;
1514  closest = h;
1515  }
1516  if ( fabs(sx) < fabs(mindx) )
1517  mindx = sx;
1518  if ( fabs(sy) < fabs(mindy) )
1519  mindy = sy;
1520 
1521  } // loop h
1522 
1523  if ( fabs(mindp) < 0.04*5 && fabs(mindr) < 9 ) {
1524  found_hits.push_back(closest);
1525  }
1526 
1527  return found_hits;
1528  } // findFttHitsNearProjectedState
1529 
1530  bool getSaveCriteriaValues() { return mSaveCriteriaValues; }
1531  std::vector<KiTrack::ICriterion *> getTwoHitCriteria() { return mTwoHitCrit; }
1532  std::vector<KiTrack::ICriterion *> getThreeHitCriteria() { return mThreeHitCrit; }
1533 
1534  TrackFitter *getTrackFitter() { return mTrackFitter; }
1535  void setEventVertex( TVector3 v, TMatrixDSym cov ){
1536  mEventVertex = v;
1537  // this is the FwdHit we will use in seeds
1538  mEventVertexHit.setXYZDetId( v.X(), v.Y(), v.Z(), kTpcId );
1539  for (size_t i=0; i < 3; i++){
1540  for (size_t j=0; j < 3; j++){
1541  mEventVertexHit._covmat(i,j) = cov(i,j);
1542  }
1543  }
1544  }
1545  TVector3 getEventVertex() { return mEventVertex; }
1546 
1547  protected:
1548  unsigned long long int nEvents;
1549 
1550  bool mDoTrackFitting = true;
1551  bool mSaveCriteriaValues = false;
1552  enum SeedSource { kFstSeed = 0, kFttSeed, kSimSeed, kSeqSeed };
1553  int mSeedSource = 1; // 0 = FST, 1 = FTT, 2 = FST+FTT simultaneous, 3 = FST+FTT sequential
1554 
1555  FwdTrackerConfig mConfig;
1556  std::string mConfigFile;
1557  size_t mTotalHitsRemoved;
1558 
1559  std::vector<GenfitTrackResult> mTrackResults;
1560 
1561  std::vector<Seed_t> mTrackSeeds; // the tracks recod from all iterations
1562  std::vector<Seed_t> mTrackSeedsThisIteration;
1563 
1564  // Metrics about the event
1565  EventStats mEventStats;
1566 
1567  // Set to the Primary vertex for the event
1568  TVector3 mEventVertex;
1569  FwdHit mEventVertexHit;
1570  genfit::GFRaveVertexFactory mGFRVertices;
1571 
1572  std::shared_ptr<FwdDataSource> mDataSource;
1573 
1574  TrackFitter *mTrackFitter = nullptr;
1575 
1576  std::vector<KiTrack::ICriterion *> mTwoHitCrit;
1577  std::vector<KiTrack::ICriterion *> mThreeHitCrit;
1578 
1579  // histograms of the raw input data
1580  TString mGeoCache;
1581 };
1582 
1583 #endif
size_t nHitsInHitMap(FwdDataSource::HitMap_t &hitmap)
Determine the total num of hits in the hitmap.
Definition: FwdTracker.h:472
void clearCriteria(std::vector< KiTrack::ICriterion * > &crits)
Clear the loaded criteria.
Definition: FwdTracker.h:349
std::vector< KiTrack::ICriterion * > loadCriteria(string path)
Loads Criteria from XML configuration. Utility function for loading criteria from XML config...
Definition: FwdTracker.h:309
void doSeedFindingIteration(size_t iIteration, FwdDataSource::HitMap_t &hitmap)
Main tracking procedure.
Definition: FwdTracker.h:1229
void findTrackSeeds()
Perform the track finding Creates a list of track seeds from the hitmaps Retrieve them using getTrack...
Definition: FwdTracker.h:556
void addFstHits(GenfitTrackResult &gtr, size_t disk)
Adds compatible FST hits to a track.
Definition: FwdTracker.h:1419
Definition: FwdHit.h:70
void doTrackFitting(const std::vector< Seed_t > &trackSeeds)
Loop on track seeds and fit each one.
Definition: FwdTracker.h:728
void removeHits(FwdDataSource::HitMap_t &hitmap, std::vector< Seed_t > &tracks)
Remove used hits from the hit map.
Definition: FwdTracker.h:486
void setSaveCriteriaValues(bool save)
Set the Save Criteria Values object.
Definition: FwdTracker.h:280
GenfitTrackResult fitTrack(Seed_t &seed, TVector3 *momentumSeedState=nullptr)
Perform a single fit from seed points.
Definition: FwdTracker.h:675
void cleanup()
cleanup the event-wise data structures
Definition: FwdTracker.h:523
void addFttHitsMc(GenfitTrackResult &gtr)
Adds compatible FTT hits using MC info.
Definition: FwdTracker.h:1384
std::vector< float > getCriteriaValues(std::string crit_name)
Get the Criteria Values object.
Definition: FwdTracker.h:365
void doMcTrackFinding(FwdDataSource::McTrackMap_t &mcTrackMap, int seedSource)
MC track finding builds track seeds from available hits using MC association.
Definition: FwdTracker.h:983
std::vector< int > getCriteriaTrackIds(std::string crit_name)
Get the Criteria Track Ids object.
Definition: FwdTracker.h:423
void setDCA(TVector3 pv)
Set the DCA and primary vertex for the event.
Definition: FwdTracker.h:180
void clearSavedCriteriaValues()
Clear the saved values for two hit and three hit criteria.
Definition: FwdTracker.h:450
virtual void initialize(TString geoCache, bool genHistograms)
Initialize FwdTracker.
Definition: FwdTracker.h:295
Seed_t findFttHitsNearProjectedState(Seed_t &available_hits, genfit::MeasuredStateOnPlane &msp, double dx=1.5, double dy=1.5)
Finds FTT hits near projected state.
Definition: FwdTracker.h:1490
void mergeHitmaps(FwdDataSource::HitMap_t &hitmap1, FwdDataSource::HitMap_t &hitmap2)
merges the FST and FTT hitmaps into a single hitmap The FTT hits are shifted by the number of FST sec...
Definition: FwdTracker.h:510
void addFstHitsMc(GenfitTrackResult &gtr)
Adds compatible FST hits to tracks seeded with FTT.
Definition: FwdTracker.h:1314
vector< Seed_t > doSeedFindingOnHitmapSubset(size_t iIteration, FwdDataSource::HitMap_t &hitmap)
Does track finding steps on a subset of hits (phi slice)
Definition: FwdTracker.h:1079
size_t sliceHitMapInPhi(FwdDataSource::HitMap_t &inputMap, FwdDataSource::HitMap_t &outputMap, float phi_min, float phi_max)
Slices a hitmap into a phi section.
Definition: FwdTracker.h:1056
void setConfigFile(std::string cf)
Set the Config File object.
Definition: FwdTracker.h:271
std::vector< std::map< std::string, float > > getCriteriaAllValues(std::string crit_name)
Get the Criteria All Values object.
Definition: FwdTracker.h:394
Seed_t findFstHitsNearProjectedState(Seed_t &available_hits, genfit::MeasuredStateOnPlane &msp, double dphi=0.004 *20.5, double dr=2.75 *2)
Finds FST hits near projected state.
Definition: FwdTracker.h:1459
void addFttHits(GenfitTrackResult &gtr, size_t disk)
Adds compatible FTT hits to tracks seeded with FST.
Definition: FwdTracker.h:1348