StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
TrackFitter.h
1 #ifndef TRACK_FITTER_H
2 #define TRACK_FITTER_H
3 
4 #include "GenFit/ConstField.h"
5 #include "GenFit/EventDisplay.h"
6 #include "GenFit/Exception.h"
7 #include "GenFit/FieldManager.h"
8 #include "GenFit/KalmanFitStatus.h"
9 #include "GenFit/GblFitter.h"
10 
11 #include "TDatabasePDG.h"
12 #include "TGeoManager.h"
13 #include "TMath.h"
14 #include "TRandom.h"
15 #include "TRandom3.h"
16 #include "TVector3.h"
17 
18 #include <vector>
19 #include <memory>
20 
21 #include "StFwdTrackMaker/Common.h"
22 
23 #include "StFwdTrackMaker/include/Tracker/FwdHit.h"
24 #include "StFwdTrackMaker/include/Tracker/TrackFitter.h"
25 #include "StFwdTrackMaker/include/Tracker/STARField.h"
26 #include "StFwdTrackMaker/include/Tracker/FwdGeomUtils.h"
27 
28 #include "StarGenerator/UTIL/StarRandom.h"
29 #include "FitterUtils.h"
30 
31 /* Class for interfacing with GenFit for fitting tracks
32  *
33  */
34 class TrackFitter {
35 
36 // Accessors and options
37  public:
38  std::shared_ptr<genfit::Track> getTrack() { return mFitTrack; }
39 
40  public:
41 
48  TrackFitter(FwdTrackerConfig _mConfig, TString geoCache) : mConfig(_mConfig), mGeoCache(geoCache), mFitTrack(nullptr) {}
49 
58  void setup() {
59 
60  // the geometry manager that GenFit will use
61  TGeoManager * gMan = nullptr;
62 
63  // Setup the Geometry used by GENFIT
64  LOG_INFO << "StFwdTrackMaker is loading the geometry cache: " << mConfig.get<string>("Geometry", mGeoCache.Data()).c_str() << endm;
65  TGeoManager::Import(mConfig.get<string>("Geometry", mGeoCache.Data()).c_str());
66  gMan = gGeoManager;
67  // Set up the material interface and set material effects on/off from the config
68  genfit::MaterialEffects::getInstance()->init(new genfit::TGeoMaterialInterface());
69 
70  // Set Material Stepper debug level
71  genfit::MaterialEffects::getInstance()->setDebugLvl( mConfig.get<int>("TrackFitter.MaterialEffects:DebugLvl", 0) );
72 
73  genfit::MaterialEffects::getInstance()->setEnergyLossBetheBloch( mConfig.get<int>("TrackFitter.MaterialEffects.EnergyLossBetheBloch", true) );
74  genfit::MaterialEffects::getInstance()->setNoiseBetheBloch( mConfig.get<int>("TrackFitter.MaterialEffects.NoiseBetheBloch", true) );
75  genfit::MaterialEffects::getInstance()->setNoiseCoulomb( mConfig.get<int>("TrackFitter.MaterialEffects.NoiseCoulomb", true) );
76  genfit::MaterialEffects::getInstance()->setEnergyLossBrems( mConfig.get<int>("TrackFitter.MaterialEffects.EnergyLossBrems", true) );
77  genfit::MaterialEffects::getInstance()->setNoiseBrems( mConfig.get<int>("TrackFitter.MaterialEffects.NoiseBrems", true) );
78  genfit::MaterialEffects::getInstance()->ignoreBoundariesBetweenEqualMaterials( mConfig.get<int>("TrackFitter.MaterialEffects.ignoreBoundariesBetweenEqualMaterials", true) );
79 
80  // do this last to override
81  genfit::MaterialEffects::getInstance()->setNoEffects( !mConfig.get<bool>("TrackFitter:MaterialEffects", false)); // negated, true means defaul is all effects on (noEffects off)
82  if (!mConfig.get<bool>("TrackFitter:MaterialEffects", false)){
83  LOG_INFO << "Turning OFF GenFit Material Effects in stepper" << endm;
84  }
85 
86  // Determine which Magnetic field to use
87  // Either constant field or real field from StarFieldAdaptor
88  if (mConfig.get<bool>("TrackFitter:constB", false)) {
89  mBField = std::unique_ptr<genfit::AbsBField>(new genfit::ConstField(0., 0., 5.0)); // 0.5 T Bz
90  LOG_INFO << "StFwdTrackMaker: Tracking with constant magnetic field" << endm;
91  } else if (mConfig.get<bool>("TrackFitter:zeroB", false)) {
92  mBField = std::unique_ptr<genfit::AbsBField>(new genfit::ConstField(0., 0., 0.)); // ZERO FIELD
93  LOG_INFO << "StFwdTrackMaker: Tracking with ZERO magnetic field" << endm;
94  } else {
95  mBField = std::unique_ptr<genfit::AbsBField>(new StarFieldAdaptor());
96  LOG_INFO << "StFwdTrackMaker: Tracking with StarFieldAdapter" << endm;
97  }
98  // we must have one of the two available fields at this point
99  // note, the pointer is still bound to the lifetime of the TackFitter
100  genfit::FieldManager::getInstance()->init(mBField.get());
101 
102  // initialize the main mFitter using a KalmanFitter with reference tracks
103  mFitter = std::unique_ptr<genfit::AbsKalmanFitter>(new genfit::KalmanFitterRefTrack());
104 
105  // Here we load several options from the config,
106  // to customize the mFitter behavior
107  mFitter->setMaxFailedHits(mConfig.get<int>("TrackFitter.KalmanFitterRefTrack:MaxFailedHits", -1)); // default -1, no limit
108  mFitter->setDebugLvl(mConfig.get<int>("TrackFitter.KalmanFitterRefTrack:DebugLvl", 0)); // default 0, no output
109  mFitter->setMaxIterations(mConfig.get<int>("TrackFitter.KalmanFitterRefTrack:MaxIterations", 40)); // default 4 iterations
110  mFitter->setMinIterations(mConfig.get<int>("TrackFitter.KalmanFitterRefTrack:MinIterations", 0)); // default 0 iterations
111 
112  // Set the fit convergence paramters
113  mFitter->setRelChi2Change( mConfig.get<double>("TrackFitter.KalmanFitterRefTrack:RelChi2Change", 1e-3) );
114  // mFitter->setAbsChi2Change( mConfig.get<double>("TrackFitter.KalmanFitterRefTrack:AbsChi2Change", 1e-3) );
115  mFitter->setDeltaPval( mConfig.get<double>("TrackFitter.KalmanFitterRefTrack:DeltaPval", 1e-3) );
116  mFitter->setBlowUpFactor( mConfig.get<double>("TrackFitter.KalmanFitterRefTrack:BlowUpFactor", 1e3) );
117 
118  // FwdGeomUtils looks into the loaded geometry and gets detector z locations if present
119  FwdGeomUtils fwdGeoUtils( gMan );
120 
121  // these default values are the default if the detector is
122  // a) not found in the geometry
123  // b) not provided in config
124 
125  // NOTE: these defaults are needed since the geometry file might not include FST (bug being worked on separately)
126  mFSTZLocations = fwdGeoUtils.fstZ(
127  mConfig.getVector<double>("TrackFitter.Geometry:fst",
128  {140.286011, 154.286011, 168.286011 }
129  // 144.633,158.204,171.271
130  )
131  );
132 
133  if ( fwdGeoUtils.fstZ( 0 ) < 1.0 ) { // returns 0.0 on failure
134  LOG_WARN << "Using FST z-locations from config or defautl, may not match hits" << endm;
135  }
136 
137  const double dzInnerFst = 1.715 + 0.04; // cm relative to "center" of disk + residual...
138  const double dzOuterFst = 0.240 + 0.04; // cm relative to "center" of disk
139 
140  // Now add the Si detector planes at the desired location
141  std::stringstream sstr;
142  sstr << "Adding FST Planes at: ";
143  string delim = "";
144  for (auto z : mFSTZLocations) {
145  mFSTPlanes.push_back(
146  genfit::SharedPlanePtr(
147  // these normals make the planes face along z-axis
148  new genfit::DetPlane(TVector3(0, 0, z), TVector3(1, 0, 0), TVector3(0, 1, 0) )
149  )
150  );
151 
152  // Inner Module FST planes
153  mFSTPlanesInner.push_back(
154  genfit::SharedPlanePtr(
155  // these normals make the planes face along z-axis
156  new genfit::DetPlane(TVector3(0, 0, z - dzInnerFst), TVector3(1, 0, 0), TVector3(0, 1, 0) )
157  )
158  );
159  mFSTPlanesInner.push_back(
160  genfit::SharedPlanePtr(
161  // these normals make the planes face along z-axis
162  new genfit::DetPlane(TVector3(0, 0, z + dzInnerFst), TVector3(1, 0, 0), TVector3(0, 1, 0) )
163  )
164  );
165  // Outer Module FST planes
166  mFSTPlanesOuter.push_back(
167  genfit::SharedPlanePtr(
168  // these normals make the planes face along z-axis
169  new genfit::DetPlane(TVector3(0, 0, z - dzOuterFst), TVector3(1, 0, 0), TVector3(0, 1, 0) )
170  )
171  );
172  mFSTPlanesOuter.push_back(
173  genfit::SharedPlanePtr(
174  // these normals make the planes face along z-axis
175  new genfit::DetPlane(TVector3(0, 0, z + dzOuterFst), TVector3(1, 0, 0), TVector3(0, 1, 0) )
176  )
177  );
178 
179  sstr << delim << z << " (-dzInner=" << z - dzInnerFst << ", +dzInner=" << z+dzInnerFst << ", -dzOuter=" << z - dzOuterFst << ", +dzOuter=" << z + dzOuterFst << ")";
180  delim = ", ";
181  }
182  LOG_DEBUG << sstr.str() << endm;
183 
184  // Now load FTT
185  // mConfig.getVector<>(...) requires a default, hence the
186  mFTTZLocations = fwdGeoUtils.fttZ(
187  mConfig.getVector<double>("TrackFitter.Geometry:ftt", {281.082,304.062,325.058,348.068})
188  );
189 
190  if ( fwdGeoUtils.fttZ( 0 ) < 1.0 ) { // returns 0.0 on failure
191  LOG_WARN << "Using FTT z-locations from config or default, may not match hits" << endm;
192  }
193 
194  if ( mFTTZLocations.size() != 4 ){
195  LOG_ERROR << "Wrong number of FTT layers, got " << mFTTZLocations.size() << " but expected 4" << endm;
196  }
197 
198  sstr.str("");
199  sstr.clear();
200  sstr << "Adding FTT Planes at: ";
201  delim = "";
202  for (auto z : mFTTZLocations) {
203  mFTTPlanes.push_back(
204  genfit::SharedPlanePtr(
205  // these normals make the planes face along z-axis
206  new genfit::DetPlane(TVector3(0, 0, z), TVector3(1, 0, 0), TVector3(0, 1, 0))
207  )
208  );
209  sstr << delim << z;
210  delim = ", ";
211  }
212  LOG_DEBUG << sstr.str() << endm;
213  }
214 
221  genfit::SharedPlanePtr getFstPlane( FwdHit * h ){
222 
223  size_t planeId = h->getSector();
224 
225  const TVector3 hitXYZ( h->getX(), h->getY(), h->getZ() );
226 
227  double phi = hitXYZ.Phi();
228  if ( phi < 0 ) phi = TMath::Pi() * 2 + phi;
229  const double phi_slice = phi / (TMath::Pi() / 6.0); // 2pi/12
230  const int phi_index = ((int)phi_slice);
231  const double r =sqrt( pow(hitXYZ.x(), 2) + pow(hitXYZ.y(), 2) );
232 
233  const size_t idx = phi_index % 2;
234  auto planeCorr = mFSTPlanesInner[planeId*2 + idx];
235  if ( r > 16 ){
236  planeCorr = mFSTPlanesOuter[planeId*2 + idx];
237  }
238  double cdz = (h->getZ() - planeCorr->getO().Z());
239 
240  if ( cdz > 0.010 ) {
241  LOG_WARN << "FST Z =" << h->getZ() << " vs CORR Plane Z = " << planeCorr->getO().Z() << " DIFF: " << cdz << " phi_slice = " << phi_slice << ", phi_index = " << phi_index << " R=" << hitXYZ.Pt() << " idx=" << idx << endm;
242  }
243 
244  return planeCorr;
245  } // GetFST PLANE
246 
253  TMatrixDSym CovMatPlane(KiTrack::IHit *h){
254  TMatrixDSym cm(2);
255  cm(0, 0) = static_cast<FwdHit*>(h)->_covmat(0, 0);
256  cm(1, 1) = static_cast<FwdHit*>(h)->_covmat(1, 1);
257  cm(0, 1) = static_cast<FwdHit*>(h)->_covmat(0, 1);
258  return cm;
259  }
260 
261 
269  genfit::MeasuredStateOnPlane projectToFst(size_t fstPlane, std::shared_ptr<genfit::Track> fitTrack) {
270  if (fstPlane > 2) {
271  genfit::MeasuredStateOnPlane nil;
272  return nil;
273  }
274 
275  auto detFst = mFSTPlanes[fstPlane];
276  // TODO: Why use 1 here?
277  genfit::MeasuredStateOnPlane tst = fitTrack->getFittedState(1);
278  // NOTE: this returns the track length if needed
279  fitTrack->getCardinalRep()->extrapolateToPlane(tst, detFst);
280 
281  return tst;
282  }
283 
291  genfit::MeasuredStateOnPlane projectToFtt(size_t iFttPlane, std::shared_ptr<genfit::Track> fitTrack) {
292  if (iFttPlane > 3) {
293  genfit::MeasuredStateOnPlane nil;
294  return nil;
295  }
296  auto fttPlane = mFTTPlanes[iFttPlane];
297  // TODO: why use 1 here?
298  genfit::MeasuredStateOnPlane tst = fitTrack->getFittedState(1);
299  // NOTE: this returns the track length if needed
300  fitTrack->getCardinalRep()->extrapolateToPlane(tst, fttPlane);
301  return tst;
302  }
303 
311  bool setupTrack(Seed_t trackSeed ) {
312 
313  // setup the track fit seed parameters
314  GenericFitSeeder gfs;
315  int seedQ = 1;
316  TVector3 seedPos(0, 0, 0);
317  TVector3 seedMom(0, 0, 10); // this default seed actually works better than a half-bad guess
318  gfs.makeSeed( trackSeed, seedPos, seedMom, seedQ );
319  LOG_DEBUG << "Setting track fit seed position = " << TString::Format( "(%f, %f, %f)", seedPos.X(), seedPos.Y(), seedPos.Z() ) << endm;
320  LOG_DEBUG << "Setting track fit seed momentum = " << TString::Format( "(%f, %f, %f)", seedMom.X(), seedMom.Y(), seedMom.Z() ) << endm;
321 
322  LOG_DEBUG << "Setting track fit seed charge = " << seedQ << endm;
323 
324  if ( seedQ == 0 ) {
325  LOG_ERROR << "Seed charge is zero, skipping track -> usually means collinear points" << endm;
326  return false;
327  }
328 
329  // create the track representations
330  // Note that multiple track reps differing only by charge results in a silent failure of GenFit
331  auto theTrackRep = new genfit::RKTrackRep(mPdgMuon * -1 * seedQ); // bc pos PDG codes are for neg particles
332 
333  // Create the track
334  mFitTrack = std::make_shared<genfit::Track>(theTrackRep, seedPos, seedMom);
335  // now add the points to the track
336 
337  int hitId(0); // hit ID
338  size_t planeId(0); // detector plane ID
339 
340  // initialize the hit coords on plane
341  TVectorD hitCoords(2);
342  hitCoords[0] = 0;
343  hitCoords[1] = 0;
344 
345  /******************************************************************************************************************
346  * loop over the hits, add them to the track
347  ******************************************************************************************************************/
348  // use these to enforce our sorting parameters
349  size_t idxFst = 0; // index of the FST hit
350  size_t idxFtt = 0; // index of the FTT hit
351  for (auto h : trackSeed) {
352  auto fh = dynamic_cast<FwdHit*>(h);
353  hitCoords[0] = h->getX();
354  hitCoords[1] = h->getY();
355 
356  /******************************************************************************************************************
357  * If the Primary vertex is included
358  ******************************************************************************************************************/
359  if ( true ) {
360  LOG_DEBUG << "Treating hit as a spacepoint" << endm;
361  if ( fh->isPV() ){
362  LOG_DEBUG << "Including primary vertex in fit" << endm;
363  }
364  TVectorD pv(3);
365  pv[0] = h->getX();
366  pv[1] = h->getY();
367  pv[2] = h->getZ();
368  LOG_INFO << "x = " << pv[0] << "+/- " << fh->_covmat(0,0) << ", y = " << pv[1] << " +/- " << fh->_covmat(1,1) << ", z = " << pv[2] << " +/- " << fh->_covmat(2,2) << endm;
369  auto tp = new genfit::TrackPoint();
370  genfit::SpacepointMeasurement *measurement = new genfit::SpacepointMeasurement(pv, fh->_covmat, fh->_detid, ++hitId, tp);
371  tp->addRawMeasurement(measurement);
372  tp->setTrack(mFitTrack.get());
373  if ( fh->isPV() ){
374  tp->setSortingParameter(0);
375  }
376  if ( fh->isFtt() ){
377  tp->setSortingParameter(4 + idxFtt);
378  idxFtt++;
379  }
380  if ( fh->isFst() ){
381  tp->setSortingParameter(1 + idxFst);
382  idxFst++;
383  }
384 
385  mFitTrack->insertPoint( tp );
386  continue;
387  }
388 
389  genfit::PlanarMeasurement *measurement = new genfit::PlanarMeasurement(hitCoords, CovMatPlane(h), fh->_detid, ++hitId, nullptr);
390 
391  planeId = h->getSector();
392  genfit::SharedPlanePtr plane;
393  if ( fh->isFtt() ){
394  planeId = fh->_vid - 9;
395  }
396  LOG_INFO << "planeId = " << planeId << ", sector " << h->getSector() << ", vid = " << fh->_vid << endm;
397  if (fh->isFtt() && mFTTPlanes.size() <= planeId) {
398  LOG_ERROR << "invalid VolumId -> out of bounds DetPlane, vid = " << dynamic_cast<FwdHit*>(h)->_vid << " vs. planeId = " << planeId << endm;
399  delete measurement;
400  continue;
401  }
402 
403  if (fh->isFtt())
404  plane = mFTTPlanes[planeId];
405  else if (fh->isFst())
406  plane = getFstPlane( fh );
407 
408  measurement->setPlane(plane, planeId);
409 
410  mFitTrack->insertPoint(new genfit::TrackPoint(measurement, mFitTrack.get()));
411  LOG_INFO << "\tsetupTrack: Hit at Z = " << h->getZ() << " with plane at Z = " << plane->getO().Z() << endm;
412 
413  if (abs(h->getZ() - plane->getO().Z()) > 0.05) {
414  LOG_WARN << "Z Mismatch h->z = " << h->getZ() << ", plane->z = "<< plane->getO().Z() <<", diff = " << abs(h->getZ() - plane->getO().Z()) << endm;
415  }
416  } // loop on trackSeed
417  return true;
418  } // setupTrack
419 
423  void performFit( std::shared_ptr<genfit::Track> t ){
424  /******************************************************************************************************************
425  * Do the fit
426  ******************************************************************************************************************/
427  try {
428 
429  // prepare the track for fitting
430  // int nFailedPoints = 0;
431  // bool changed = false;
432  // changed = dynamic_cast<genfit::KalmanFitterRefTrack*>( mFitter.get() )->prepareTrack( mFitTrack.get(), mFitTrack->getCardinalRep(), false, nFailedPoints);
433  // LOG_DEBUG << "Track prepared for fit with " << nFailedPoints << " failed points, changed? = " << changed << endm;
434 
435  // check the track for consistency
436  mFitTrack->checkConsistency();
437  // do the fit
438  mFitter->processTrack(t.get());
439 
440  // check the track for consistency
441  t->checkConsistency();
442 
443  // find track rep with smallest chi2
444  t->determineCardinalRep();
445  // update the seed
446  // t->udpateSeed();
447 
448  auto status = t->getFitStatus();
449  LOG_INFO << "Fit status: " << status->isFitConverged() << endm;
450  LOG_INFO << "-Fit pvalue: " << status->getPVal() << endm;
451  LOG_INFO << "-Fit Chi2: " << status->getChi2() << endm;
452 
453  if ( status->isFitConverged() ){
454 
455  auto cr = t->getCardinalRep();
456  auto p = cr->getMom( t->getFittedState( 0, cr ));
457  int rcQ = status->getCharge();
458  LOG_INFO << "Fit momentum: " << p.X() << ", " << p.Y() << ", " << p.Z() << endm;
459  LOG_INFO << "\tFit Pt: " << p.Pt() << ", eta: " << p.Eta() << ", phi: " << p.Phi() << endm;
460  // LOG_INFO << "\tMc Pt: " << mcMom.Pt() << ", eta: " << mcMom.Eta() << ", phi: " << mcMom.Phi() << endm;
461  }
462 
463 
464  } catch (genfit::Exception &e) {
465  LOG_ERROR << "Exception on fit update" << e.what() << endm;
466  }
467  LOG_INFO << "Track fit update complete!" << endm;
468  }
469 
478  long long fitTrack(Seed_t trackSeed, TVector3 *seedMomentum = 0) {
479  long long itStart = FwdTrackerUtils::nowNanoSecond();
480  LOG_DEBUG << "Fitting track with " << trackSeed.size() << " FWD Measurements" << endm;
481 
482  /******************************************************************************************************************
483  * First sort the seed, bc GENFIT seemingly cannot handle out of order points
484  ******************************************************************************************************************/
485  std::sort(trackSeed.begin(), trackSeed.end(),
486  [](KiTrack::IHit *a, KiTrack::IHit *b)
487  { return a->getZ() < b->getZ(); }
488  );
489 
490  /******************************************************************************************************************
491  * Setup the track fit seed parameters and objects
492  ******************************************************************************************************************/
493  bool valid = setupTrack(trackSeed);
494  if ( !valid ){
495  LOG_ERROR << "Failed to setup track for fit" << endm;
496  return -1;
497  }
498  LOG_DEBUG << "Ready to fit with " << mFitTrack->getNumPoints() << " track points" << endm;
499 
500  /******************************************************************************************************************
501  * Do the fit
502  ******************************************************************************************************************/
503  performFit( mFitTrack );
504  long long duration = (FwdTrackerUtils::nowNanoSecond() - itStart) * 1e-6; // milliseconds
505  return duration;
506  } // fitTrack
507 
508  // Store the planes for FTT and FST
509  vector<genfit::SharedPlanePtr> mFTTPlanes;
510  vector<genfit::SharedPlanePtr> mFSTPlanes;
511  vector<genfit::SharedPlanePtr> mFSTPlanesInner;
512  vector<genfit::SharedPlanePtr> mFSTPlanesOuter;
513 
514  protected:
515  std::unique_ptr<genfit::AbsBField> mBField;
516 
517  FwdTrackerConfig mConfig; // main config object
518  TString mGeoCache;
519 
520  // Main GenFit fitter instance
521  std::unique_ptr<genfit::AbsKalmanFitter> mFitter = nullptr;
522 
523  // PDG codes for the default plc type for fits
524  const int mPdgPiPlus = 211;
525  const int mPdgPiMinus = -211;
526  const int mPdgPositron = 11;
527  const int mPdgElectron = -11;
528  const int mPdgMuon = 13;
529  const int mPdgAntiMuon = -13;
530 
531 
532  // det z locations loaded from geom or config
533  vector<double> mFSTZLocations, mFTTZLocations;
534 
535  // GenFit state - resused
536  std::shared_ptr<genfit::Track> mFitTrack;
537 };
538 
539 #endif
void performFit(std::shared_ptr< genfit::Track > t)
performs the fit on a track
Definition: TrackFitter.h:423
genfit::SharedPlanePtr getFstPlane(FwdHit *h)
Get the Fst Plane object for a given hit.
Definition: TrackFitter.h:221
Definition: FwdHit.h:70
genfit::MeasuredStateOnPlane projectToFtt(size_t iFttPlane, std::shared_ptr< genfit::Track > fitTrack)
Get projection to given FTT plane.
Definition: TrackFitter.h:291
TMatrixDSym CovMatPlane(KiTrack::IHit *h)
Convert the 3x3 covmat to 2x2 by dropping z.
Definition: TrackFitter.h:253
std::vector< T > getVector(std::string path, std::vector< T > dv) const
Get a Vector object from config.
genfit::MeasuredStateOnPlane projectToFst(size_t fstPlane, std::shared_ptr< genfit::Track > fitTrack)
Get projection to given FST plane.
Definition: TrackFitter.h:269
TrackFitter(FwdTrackerConfig _mConfig, TString geoCache)
Construct a new Track Fitter object.
Definition: TrackFitter.h:48
T get(std::string path, T dv) const
template function for getting any type that can be converted from string via stringstream ...
bool setupTrack(Seed_t trackSeed)
setup the track from the given seed and optional primary vertex
Definition: TrackFitter.h:311
long long fitTrack(Seed_t trackSeed, TVector3 *seedMomentum=0)
Primary track fitting routine.
Definition: TrackFitter.h:478
void setup()
Setup the tracker object Load geometry Setup Material Effects Setup the magnetic field Setup the fitt...
Definition: TrackFitter.h:58