StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFwdTrackMaker.h
1 #ifndef ST_FWD_TRACK_MAKER_H
2 #define ST_FWD_TRACK_MAKER_H
3 
4 #include "StChain/StMaker.h"
5 
6 #ifndef __CINT__
7 #include "GenFit/Track.h"
8 #include "StFwdTrackMaker/include/Tracker/FwdHit.h"
9 #endif
10 
11 #include "FwdTrackerConfig.h"
12 #include "TVector3.h"
13 #include "TMatrix.h"
14 
15 namespace KiTrack {
16 class IHit;
17 };
18 
19 namespace genfit {
20  class Track;
21  class GFRaveVertex;
22 }
23 
24 class ForwardTracker;
25 class ForwardTrackMaker;
26 class FwdDataSource;
27 class FwdHit;
28 class StarFieldAdaptor;
29 
30 class StGlobalTrack;
31 class StRnDHitCollection;
32 class StTrack;
34 class SiRasterizer;
35 class McTrack;
36 
37 class StFcsDb;
38 
39 // ROOT includes
40 #include "TNtuple.h"
41 #include "TTree.h"
42 // STL includes
43 #include <vector>
44 #include <memory>
45 class StFwdTrack;
46 class GenfitTrackResult;
47 
48 
49 class StFwdTrackMaker : public StMaker {
50 
51  ClassDef(StFwdTrackMaker, 0);
52 
53  public:
55  ~StFwdTrackMaker(){/* nada */};
56 
57  int Init();
58  int Finish();
59  int Make();
60  void Clear(const Option_t *opts = "");
61 
62  enum { kInnerGeometry,
63  kOuterGeometry };
64 
65  void SetConfigFile(std::string n) {
66  mConfigFile = n;
67  LoadConfiguration();
68  }
69  void LoadConfiguration();
70  void SetVisualize( bool _viz ) { mVisualize = _viz; }
71 
72  vector<StFwdTrack*> mFwdTracks;
73 
74  vector<FwdHit> &GetFttHits() { return mFwdHitsFtt; }
75  vector<FwdHit> &GetFstHits() { return mFwdHitsFst; }
76 
77  #ifndef __CINT__
78  // Get the FwdTracker object
79  std::shared_ptr<ForwardTracker> GetForwardTracker() { return mForwardTracker; }
80  const std::vector<Seed_t> &getTrackSeeds() const;
81  const std::vector<GenfitTrackResult> &getFitResults() const;
82  #endif
83  TVector3 GetEventPrimaryVertex();
84 
85  private:
86  protected:
87 
88  StFcsDb* mFcsDb = 0; // Pointer to fcs db object
89 
90  // Event Filters
91  float mEventFilterMinTofMult = 2;
92  bool mEventFilterRequireEventVertex = false;
93  bool mEventFilterRequireVpdVertex = true;
94  float mEventFilterMinVpdZ = -99;
95  float mEventFilterMaxVpdZ = 99;
96 
97  // for Wavefront OBJ export
98  size_t eventIndex = 0; // counts up for processed events
99  size_t mEventNum = 0; // global event num (index)
100  TVector3 mEventVertex; // primary vertex used in fwd tracking this event
101 
102  std::string mConfigFile;
103 
104  std::map<std::string, TH1 *> mHistograms;
105 
106  bool mVisualize = false; // if true,write out a Wavefront OBJ to visualize the event in 3D
107  vector<TVector3> mFttHits;
108  vector<TVector3> mFstHits;
109  vector<TVector3> mFcsClusters;
110  vector<float> mFcsClusterEnergy;
111  vector<TVector3> mFcsPreHits;
112 
113  std::vector< genfit::GFRaveVertex * > mRaveVertices;
114  vector<float> mFttZFromGeom, mFstZFromGeom;
115 
116  void ProcessFwdTracks();
117  void FillEvent();
118  void FillTrackDeltas();
119  bool SkipEvent();
120 
121  StFwdTrack * makeStFwdTrack( GenfitTrackResult &gtr, size_t indexTrack );
122 
123  // I could not get the library generation to succeed with these.
124  // so I have removed them
125  #ifndef __CINT__
126  TMatrixDSym mEventVertexCov; // covariance matrix for the primary vertex
127  enum FwdVertexSource { kFwdVertexSourceUnknown, kFwdVertexSourceNone, kFwdVertexSourceTpc, kFwdVertexSourceMc, kFwdVertexSourceVpd }; // unknown means we havent looked yet
128  FwdVertexSource mFwdVertexSource = StFwdTrackMaker::kFwdVertexSourceUnknown;
129  vector<FwdHit> mFwdHitsFtt;
130  vector<FwdHit> mFwdHitsFst;
131  std::shared_ptr<SiRasterizer> mSiRasterizer;
132  FwdTrackerConfig mFwdConfig;
133  std::shared_ptr<ForwardTracker> mForwardTracker;
134  std::shared_ptr<FwdDataSource> mForwardData;
135  size_t loadMcTracks( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap );
136  void loadFcs();
137  void loadFttHits( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap, int count = 0 );
138  void loadFttHitsFromStEvent( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap, int count = 0 );
139  void loadFttHitsFromGEANT( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap, int count = 0 );
140 
141  int loadFstHits( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap );
142  int loadFstHitsFromMuDst( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap );
143  int loadFstHitsFromGEANT( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap );
144  int loadFstHitsFromStEvent( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap );
145  int loadFstHitsFromStRnDHits( std::map<int, std::shared_ptr<McTrack>> &mcTrackMap, std::map<int, std::vector<KiTrack::IHit *>> &hitMap );
146  #endif
147 
148 
150  void FitVertex();
151 
152  static std::string defaultConfig;
153  bool configLoaded = false;
154  TString mGeoCache;
155 
156  // Helper functions for modifying configuration
157  public:
161  void setOutputFilename( std::string fn ) { mFwdConfig.set( "Output:url", fn ); }
166  void setFttHitSource( std::string source ) { mFwdConfig.set( "Source:ftt", source ); }
167 
171  void setUseFstRasteredGeantHits( bool use = true ){ mFwdConfig.set<bool>( "SiRasterizer:active", use ); }
176  void setFstRasterR( double r = 3.0 /*cm*/ ){ mFwdConfig.set<double>( "SiRasterizer:r", r ); }
181  void setFstRasterPhi( double phi = 0.00409 /*2*pi/(12*128)*/ ){ mFwdConfig.set<double>( "SiRasterizer:phi", phi ); }
182 
183  //Track Finding
187  void setSeedFindingWithFstFttSequential() { mFwdConfig.set( "TrackFinder:source", "seq" ); }
191  void setSeedFindingWithFstFttSimultaneous() { mFwdConfig.set( "TrackFinder:source", "sim" ); }
195  void setSeedFindingWithFtt() { mFwdConfig.set( "TrackFinder:source", "ftt" ); }
199  void setSeedFindingWithFst() { mFwdConfig.set( "TrackFinder:source", "fst" ); }
203  void setSeedFindingNumInterations( int n = 1 ) { mFwdConfig.set<int>("TrackFinder:nIterations", n); }
207  void setSeedFindingNumPhiSlices( int n = 8 ) { mFwdConfig.set<int>("TrackFinder.Iteration:nPhiSlices", n); }
211  void setSeedFindingConnectorDistance( int d = 1 ) { mFwdConfig.set<int>( "TrackFinder.Connector:distance", d ); }
216  void setSeedFindingUseSubsetNN( bool use = true ) { mFwdConfig.set<bool>( "TrackFinder.SubsetNN:active", use ); }
220  void setSeedFindingMinHitsOnTrack( int n = 3 ) { mFwdConfig.set<int>( "TrackFinder.SubsetNN:min-hits-on-track", n ); }
225  void setSeedFindingUseHitRemover( bool use = true ) { mFwdConfig.set<bool>( "TrackFinder.HitRemover:active", use ); }
230  void setUseTruthSeedFinding( bool use = true ) { mFwdConfig.set<bool>( "TrackFinder:active", !use ); }
231 
232  // Track Fitting
236  void setTrackFittingOff() { mFwdConfig.set( "TrackFitter:active", "false" ); }
240  void setFittingMaterialEffects( bool mat = true) { mFwdConfig.set<bool>( "TrackFitter:materialEffects", mat ); }
244  void setPrimaryVertexSigmaXY( double sXY ) { mFwdConfig.set<double>( "TrackFitter.Vertex:sigmaXY", sXY ); }
248  void setPrimaryVertexSigmaZ( double sZ ) { mFwdConfig.set<double>( "TrackFitter.Vertex:sigmaZ", sZ ); }
249  // TODO: add options for beamline constraint
253  void setZeroB( bool zeroB = true ) { mFwdConfig.set<bool>( "TrackFitter:zeroB", zeroB ); }
257  void setConstB( bool constB = true ) { mFwdConfig.set<bool>( "TrackFitter:constB", constB ); }
261  void setUseMcSeedForFit( bool mcSeed = true ) { mFwdConfig.set<bool>( "TrackFitter:mcSeed", mcSeed ); }
262 
269  void setTrackRefit( bool refit = true) { mFwdConfig.set<bool>( "TrackFitter:refit", refit ); }
270 
274  void setMaxFailedHitsInFit( int n = -1 /*no lim*/ ) {mFwdConfig.set<int>("TrackFitter.KalmanFitterRefTrack:MaxFailedHits", n);}
278  void setFitDebugLvl( int level = 0 /*0=no output*/ ) {mFwdConfig.set<int>("TrackFitter.KalmanFitterRefTrack:DebugLvl", level); }
282  void setFitMaxIterations( int n=4 ) {mFwdConfig.set<int>("TrackFitter.KalmanFitterRefTrack:MaxIterations", n); }
286  void setFitMinIterations( int n = 1) {mFwdConfig.set<int>("TrackFitter.KalmanFitterRefTrack:MinIterations", n); }
287 
291  void setSmearMcPrimaryVertex( bool pvs = true ) { mFwdConfig.set<bool>( "TrackFitter.Vertex:smearMcVertex", pvs ); }
292 
297  void setGeoCache( TString gc ) { mGeoCache = gc; }
298 
305  void setConfigKeyValue( std::string k, std::string v ){
306  mFwdConfig.set( k, v );
307  }
308 
314  void setCrit2( std::string name, double min, double max ){
315  for ( auto p : mFwdConfig.childrenOf( "TrackFinder.Iteration.SegmentBuilder" ) ){
316  auto nName = mFwdConfig.get<std::string>( p + ":name", "DNE" );
317  if (nName == name) {
318  LOG_DEBUG << "Setting Crit2=" << nName << " (min=" << min << ", max=" << max << ")" << endm;
319  mFwdConfig.set<double>(p + ":min", min );
320  mFwdConfig.set<double>(p + ":max", max );
321  return;
322  }
323  } // loop on existing crit2
324  // if we got here then the crit did not exist
325 
326  }
327 
333  void setCrit3( std::string name, double min, double max ){
334  for ( auto p : mFwdConfig.childrenOf( "TrackFinder.Iteration.ThreeHitSegments" ) ){
335  auto nName = mFwdConfig.get<std::string>( p + ":name", "DNE" );
336  if (nName == name) {
337  LOG_DEBUG << "Setting Crit3=" << nName << " (min=" << min << ", max=" << max << ")" << endm;
338  mFwdConfig.set<double>(p + ":min", min );
339  mFwdConfig.set<double>(p + ":max", max );
340  return;
341  }
342  } // loop on existing crit3
343  // if we got here then the crit did not exist
344  }
345 
346 };
347 
348 #endif
void setCrit3(std::string name, double min, double max)
Sets a criteria value in the config for 3-hit criteria.
void setSmearMcPrimaryVertex(bool pvs=true)
Enables smearing of the MC Primary Vertex according to sigmaXY,Z.
int loadFstHits(std::map< int, std::shared_ptr< McTrack >> &mcTrackMap, std::map< int, std::vector< KiTrack::IHit * >> &hitMap)
Loads FST hits from various sources into the hitmap and McTrackMap (if availabale) ...
void setZeroB(bool zeroB=true)
Set B-field to zero (for zero field running)
void FitVertex()
Fit the primary vertex using FWD tracks.
void setPrimaryVertexSigmaZ(double sZ)
Set the resolution for the Primary Vertex in Z sZ : sigma in Z (cm)
void setSeedFindingUseHitRemover(bool use=true)
Enable or disable the HitRemover.
void setFittingMaterialEffects(bool mat=true)
Enable / disable material effects Material effects in kalman filter.
void setConfigKeyValue(std::string k, std::string v)
Set a generic Key Value in the Config object.
void setCrit2(std::string name, double min, double max)
Sets a criteria value in the config for 2-hit criteria.
Definition: FwdHit.h:70
void setFttHitSource(std::string source)
Set the data source for FTT hits.
std::vector< std::string > childrenOf(std::string path) const
list the paths of children nodes for a given node
void setUseMcSeedForFit(bool mcSeed=true)
Force the use of McSeed for fit.
StFwdTrack * makeStFwdTrack(GenfitTrackResult &gtr, size_t indexTrack)
void setMaxFailedHitsInFit(int n=-1)
Sets the maximum number of hits that can be considered failed before the entire track fit fails...
void setSeedFindingWithFstFttSimultaneous()
Use FST and Ftt hits (simultaneously) in the Seed Finding.
void setSeedFindingConnectorDistance(int d=1)
Set the connector distance for track finding.
void setUseFstRasteredGeantHits(bool use=true)
Enable or disable the Fst Rasterizer.
void setSeedFindingUseSubsetNN(bool use=true)
Enable or disable the SubsetNN.
void setSeedFindingWithFstFttSequential()
Use FST and Ftt hits (sequentially) in the Seed Finding - then merge tracks.
void setConstB(bool constB=true)
Set B-field to constant (even outside of TPC)
void setTrackFittingOff()
Turn off track fitting Useful if you want to speed up the run but dont need fitting (testing seed fin...
Definition: FwdHit.h:37
void setSeedFindingWithFtt()
Use Ftt hits in the Seed Finding.
void set(std::string path, T v)
Writes a value of type T to the map Uses convertTo&lt;T&gt; to convert type T to a string rep...
void setSeedFindingNumInterations(int n=1)
Set the number of track finding iterations.
void setFstRasterPhi(double phi=0.00409)
Set the resolution in phi for rasterizing FST hits (from fast sim) Only used when the Rasterizer is e...
void setFitMaxIterations(int n=4)
Sets Max fit iterations before failing.
void setFstRasterR(double r=3.0)
Set the resolution in R for rasterizing FST hits (from fast sim) Only used when the Rasterizer is ena...
void setSeedFindingMinHitsOnTrack(int n=3)
Enable or disable the SubsetNN.
void setUseTruthSeedFinding(bool use=true)
Enable or disable the Truth Seed finding.
void setOutputFilename(std::string fn)
Set the filename for output ROOT file.
T get(std::string path, T dv) const
template function for getting any type that can be converted from string via stringstream ...
void setFitDebugLvl(int level=0)
Sets Fitter debug level.
void setTrackRefit(bool refit=true)
Sets the tracking to refit This adds compatible hits from whichever detector was NOT used in seed fin...
size_t loadMcTracks(std::map< int, std::shared_ptr< McTrack >> &mcTrackMap)
void setGeoCache(TString gc)
Sets geometry cache filename.
void setSeedFindingNumPhiSlices(int n=8)
Set the number of phi slices to split the track iterations into.
void setFitMinIterations(int n=1)
Sets Min fit iterations before converging.
void setPrimaryVertexSigmaXY(double sXY)
Set the resolution for the Primary Vertex in XY sXY : sigma in XY (cm)
C++ STL includes.
Definition: AgUStep.h:47
void setSeedFindingWithFst()
Use Fst hits in the Seed Finding.