FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndRecoKalmanFit Class Reference

#include <PndRecoKalmanFit.h>

Inheritance diagram for PndRecoKalmanFit:

Public Member Functions

 PndRecoKalmanFit ()
 
 ~PndRecoKalmanFit ()
 
void SetGeane (Bool_t opt=kTRUE)
 
void SetPropagateToIP (Bool_t opt=kTRUE)
 
void SetPropagateDistance (Float_t opt=-1.f)
 
void SetPerpPlane (Bool_t opt=kTRUE)
 
void SetNumIterations (Int_t num)
 
void SetTrackRep (Int_t num)
 
void SetVerbose (Int_t verb)
 
void SetMvdBranchName (const TString &name)
 
void SetCentralTrackerBranchName (const TString &name)
 
Bool_t Init ()
 
PndTrackFit (PndTrack *tBefore, Int_t PDG)
 
GFRecoHitFactoryGetRecoHitFactory ()
 

Private Member Functions

 ClassDef (PndRecoKalmanFit, 0)
 Verbose level. More...
 

Private Attributes

GFRecoHitFactoryfTheRecoHitFactory
 
GFKalman fGenFitter
 
FairGeanePro * fPro
 
TString fMvdBranchName
 Geane Propagator. More...
 
TString fCentralTrackerBranchName
 Name of the TCA for MVD. More...
 
Bool_t fUseGeane
 Name of the TCA for central tracker. More...
 
Bool_t fPropagateToIP
 Flag to use Geane. More...
 
Float_t fPropagateDistance
 Flag to propagate to the interaction point. More...
 
Bool_t fPerpPlane
 Distance in [cm] to back-propagate the parameters, negative number means no backpropagation. More...
 
Int_t fNumIt
 Flag to use as initial plane the one perpendicular to the track. More...
 
Short_t fTrackRep
 Number of iterations. More...
 
Int_t fVerbose
 (0) GeaneTrackRep, 1 RKTrackRep More...
 

Detailed Description

Definition at line 31 of file PndRecoKalmanFit.h.

Constructor & Destructor Documentation

PndRecoKalmanFit::PndRecoKalmanFit ( )

Definition at line 70 of file PndRecoKalmanFit.cxx.

References PndGeoHandling::Instance().

70  : TNamed("Genfit", "Fit Tracks"),
72  fUseGeane(kTRUE), fPropagateToIP(kTRUE), fPropagateDistance(-1.f), fPerpPlane(kFALSE), fNumIt(1), fTrackRep(0), fVerbose(0)
73 {
75 }
TString fCentralTrackerBranchName
Name of the TCA for MVD.
Short_t fTrackRep
Number of iterations.
Bool_t fPropagateToIP
Flag to use Geane.
TString fMvdBranchName
Geane Propagator.
Bool_t fPerpPlane
Distance in [cm] to back-propagate the parameters, negative number means no backpropagation.
TFile * f
Definition: bump_analys.C:12
Int_t fVerbose
(0) GeaneTrackRep, 1 RKTrackRep
static PndGeoHandling * Instance()
Float_t fPropagateDistance
Flag to propagate to the interaction point.
Int_t fNumIt
Flag to use as initial plane the one perpendicular to the track.
Bool_t fUseGeane
Name of the TCA for central tracker.
PndRecoKalmanFit::~PndRecoKalmanFit ( )

Definition at line 205 of file PndRecoKalmanFit.cxx.

References fPro.

205  {
206  delete(fPro);
207 }
FairGeanePro * fPro

Member Function Documentation

PndRecoKalmanFit::ClassDef ( PndRecoKalmanFit  ,
 
)
private

Verbose level.

PndTrack * PndRecoKalmanFit::Fit ( PndTrack tBefore,
Int_t  PDG 
)

Definition at line 209 of file PndRecoKalmanFit.cxx.

References GFTrack::addHitVector(), Bool_t, GFRecoHitFactory::createMany(), f, fabs(), fGenFitter, fPerpPlane, fPro, fPropagateDistance, fPropagateToIP, fTheRecoHitFactory, fTrackRep, fVerbose, GenfitTrack2PndTrack(), GFTrack::getCand(), GFTrack::getNumHits(), PndTrack::GetParamFirst(), PndTrack::GetTrackCand(), par, PndTrackCand2GenfitTrackCand(), GFKalman::processTrack(), GFTrack::setCandidate(), PndTrack::SetFlag(), GeaneTrackRep::setPropDir(), trk, and GFException::what().

Referenced by PndRecoMultiKalmanTask::Exec(), and PndRecoKalmanTask::Exec().

210 {
211  PndTrack* tAfter = NULL;
212  if (fVerbose>0) std::cout<<"PndRecoKalmanFit::Fit"<<std::endl;
213  if (fabs(tBefore->GetParamFirst().GetPz())<1e-9)
214  {
215  tAfter = tBefore;
216  tAfter->SetFlag(-10);
217  return tAfter; // flag -10 : pz==0
218  }
219 
220  Int_t fCharge= tBefore->GetParamFirst().GetQ();
221  Int_t PDGCode= PDG;
222  TVector3 StartPos(tBefore->GetParamFirst().GetX(),tBefore->GetParamFirst().GetY(),tBefore->GetParamFirst().GetZ());
223  TVector3 StartMom(tBefore->GetParamFirst().GetPx(),tBefore->GetParamFirst().GetPy(),tBefore->GetParamFirst().GetPz());
224  TVector3 StartPosErr(tBefore->GetParamFirst().GetDX(),tBefore->GetParamFirst().GetDY(),tBefore->GetParamFirst().GetDZ());
225  TVector3 StartMomErr(tBefore->GetParamFirst().GetDPx(),tBefore->GetParamFirst().GetDPy(),tBefore->GetParamFirst().GetDPz());
226 
227  GFAbsTrackRep* rep = 0;
228 
229  if (fPropagateToIP)
230  {
231  // Calculating params at PCA to Origin
232  FairTrackParP par = tBefore->GetParamFirst();
233  Int_t ierr = 0;
234  FairTrackParH helix(&par, ierr);
235  FairGeanePro fPro0;
236  if (fVerbose==0) fPro0.SetPrintErrors(kFALSE);
237  FairTrackParH fRes;
238  fPro0.SetPoint(TVector3(0,0,0));
239  fPro0.PropagateToPCA(1, -1);
240  Bool_t rc = fPro0.Propagate(&helix, &fRes, PDGCode);
241  if (rc)
242  {
243  StartPos.SetXYZ(fRes.GetX(), fRes.GetY(), fRes.GetZ());
244  StartMom.SetXYZ(fRes.GetPx(), fRes.GetPy(), fRes.GetPz());
245  StartPosErr.SetXYZ(fRes.GetDX(), fRes.GetDY(), fRes.GetDZ());
246  StartMomErr.SetXYZ(fRes.GetDPx(), fRes.GetDPy(), fRes.GetDPz());
247  }
248  }
249  else if (fPropagateDistance>0.f)
250  {
251  // Calculating params at fPropagateDistance cm before the first hit
252  FairTrackParP par = tBefore->GetParamFirst();
253  Int_t ierr = 0;
254  FairTrackParH helix(&par, ierr);
255  FairGeanePro fPro0;
256  if (fVerbose==0) fPro0.SetPrintErrors(kFALSE);
257  FairTrackParH fRes;
258  fPro0.PropagateToLength(-fPropagateDistance);
259  Bool_t rc = fPro0.Propagate(&helix, &fRes, PDGCode);
260  if (rc)
261  {
262  StartPos.SetXYZ(fRes.GetX(), fRes.GetY(), fRes.GetZ());
263  StartMom.SetXYZ(fRes.GetPx(), fRes.GetPy(), fRes.GetPz());
264  StartPosErr.SetXYZ(fRes.GetDX(), fRes.GetDY(), fRes.GetDZ());
265  StartMomErr.SetXYZ(fRes.GetDPx(), fRes.GetDPy(), fRes.GetDPz());
266  }
267  }
268 
269  TVector3 plane_v1, plane_v2;
270  if (fPerpPlane)
271  {
272  plane_v1 = StartMom.Orthogonal();
273  plane_v2 = StartMom.Cross(plane_v1);
274  }
275  else
276  {
277  plane_v1.SetXYZ(1.,0.,0.);
278  plane_v2.SetXYZ(0.,1.,0.);
279  }
280  GFDetPlane start_pl(StartPos, plane_v1, plane_v2);
281  GFTrack* trk;
282  if (fTrackRep==0)
283  {
284  GeaneTrackRep *grep = new GeaneTrackRep(fPro,
285  start_pl,StartMom,
286  StartPosErr,StartMomErr,
287  fCharge,PDGCode);
288  grep->setPropDir(1);
289  rep = grep;
290  }
291  else if (fTrackRep==1)
292  {
293  RKTrackRep *grep = new RKTrackRep(StartPos, StartMom,
294  StartPosErr, StartMomErr,
295  PDGCode);
296  rep = grep;
297  }
298  else
299  {
300  std::cout << "*** PndRecoKalmanFit::Fit" << "\t" << "Not existing Track Representation " << fTrackRep << std::endl;
301  return NULL; // any smarted ideas?
302  }
303 
304  trk= new GFTrack(rep);
305 
306  PndTrackCand trackCand = tBefore->GetTrackCand();
307  trk->setCandidate(*PndTrackCand2GenfitTrackCand(&trackCand));
308 
309  // Load RecoHits
310  try
311  {
313  if (fVerbose>0) std::cout<<trk->getNumHits()<<" hits in track " << std::endl;
314  }
315  catch(GFException& e)
316  {
317  std::cout << "*** PndRecoKalmanFit::Fit" << "\t" << "Genfit Exception: trk->addHitVector " << e.what() << std::endl;
318  //throw e;
319  }
320  // Start Fitter
321  try
322  {
324  }
325  catch (GFException e)
326  {
327  std::cout<<"*** PndRecoKalmanFit::Fit" << "\t" << "FITTER EXCEPTION ***"<<std::endl;
328  std::cout<<e.what()<<std::endl;
329  }
330  if (fVerbose>0) std::cout<<"*** PndRecoKalmanFit::Fit" << "\t" << "SUCCESSFULL FIT!"<<std::endl;
331 
332  try
333  {
334  tAfter = (PndTrack*)GenfitTrack2PndTrack(trk);
335  }
336  catch (GFException e)
337  {
338  std::cout<<"*** PndRecoKalmanFit::Fit" << "\t" << "CONVERSION EXCEPTION ***"<<std::endl;
339  std::cout<<e.what()<<std::endl;
340  tAfter = tBefore;
341  tAfter->SetFlag(-2); // flag -2: conversion failed
342  }
343 
344  if (fVerbose>0) std::cout<<"Fitting done"<<std::endl;
345 
346  delete(trk);
347 
348  return tAfter;
349 }
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:80
unsigned int getNumHits() const
Definition: GFTrack.h:148
Detector plane genfit geometry class.
Definition: GFDetPlane.h:59
Track object for genfit. genfit algorithms work on these objects.
Definition: GFTrack.h:60
std::vector< GFAbsRecoHit * > createMany(const GFTrackCand &cand)
Creat a collection of RecoHits.
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:47
Short_t fTrackRep
Number of iterations.
Double_t par[3]
void addHitVector(std::vector< GFAbsRecoHit * > hits)
Add collection of hits.
Definition: GFTrack.h:310
FairGeanePro * fPro
Bool_t fPropagateToIP
Flag to use Geane.
void setCandidate(const GFTrackCand &cand, bool doreset=false)
set track candidate
Definition: GFTrack.cxx:148
Track Representation module based on a Runge-Kutta algorithm including a full material model...
Bool_t fPerpPlane
Distance in [cm] to back-propagate the parameters, negative number means no backpropagation.
PndTrackCand GetTrackCand()
Definition: PndTrack.h:47
TFile * f
Definition: bump_analys.C:12
void setPropDir(int d)
GFTrackCand * PndTrackCand2GenfitTrackCand(PndTrackCand *cand)
void processTrack(GFTrack *trk)
Performs fit on a GFTrack.
Definition: GFKalman.cxx:40
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
Int_t fVerbose
(0) GeaneTrackRep, 1 RKTrackRep
Float_t fPropagateDistance
Flag to propagate to the interaction point.
GFTrack * trk
Definition: checkgenfit.C:13
const GFTrackCand & getCand() const
Definition: GFTrack.h:142
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:45
PndTrack * GenfitTrack2PndTrack(const GFTrack *tr)
void SetFlag(Int_t i)
Definition: PndTrack.h:38
GFRecoHitFactory * fTheRecoHitFactory
FairTrackParP GetParamFirst()
Definition: PndTrack.h:49
GFRecoHitFactory* PndRecoKalmanFit::GetRecoHitFactory ( )
inline

Definition at line 53 of file PndRecoKalmanFit.h.

References fTheRecoHitFactory.

53 { return fTheRecoHitFactory;};
GFRecoHitFactory * fTheRecoHitFactory
Bool_t PndRecoKalmanFit::Init ( )

Definition at line 77 of file PndRecoKalmanFit.cxx.

References GFRecoHitFactory::addProducer(), fCentralTrackerBranchName, fGenFitter, PndFtsMapCreator::FillTubeArray(), PndSttMapCreator::FillTubeArray(), fMvdBranchName, fNumIt, fPro, fTheRecoHitFactory, fTrackRep, fUseGeane, fVerbose, PndGeoSttPar::GetGeometryType(), PndGeoFtsPar::GetGeometryType(), GFFieldManager::getInstance(), GFFieldManager::init(), GFException::quiet(), rtdb, and GFKalman::setNumIterations().

Referenced by PndRecoMultiKalmanTask::Init(), PndRecoKalmanTask::Init(), and PndPidCorrelator::Init().

78 {
79  //Get ROOT Manager
80  FairRootManager* ioman= FairRootManager::Instance();
81  if(ioman==0)
82  {
83  Error("PndRecoKalmanFit::Init","RootManager not instantiated!");
84  return kFALSE;
85  }
86 
87  // STT map loading
88  FairRuntimeDb* rtdb = FairRunAna::Instance()->GetRuntimeDb();
89  PndGeoSttPar *sttParameters = (PndGeoSttPar*) rtdb->getContainer("PndGeoSttPar");
90  TClonesArray *tubeArray = NULL;
91  if(sttParameters->GetGeometryType() != -1) {
92  PndSttMapCreator mapper(sttParameters);
93  tubeArray = mapper.FillTubeArray();
94  }
95  // FTS map loading
96  PndGeoFtsPar *ftsParameters = (PndGeoFtsPar*) rtdb->getContainer("PndGeoFtsPar");
97  TClonesArray *ftsTubeArray = NULL;
98  if(ftsParameters->GetGeometryType() != -1) {
99  PndFtsMapCreator ftsMapper(ftsParameters);
100  ftsTubeArray = ftsMapper.FillTubeArray();
101  }
102 
103 
104  // Build hit factory -----------------------------
106  if (fVerbose<2) GFException::quiet(true);
107 
108  TClonesArray* stripar; TClonesArray* pixelar;
109  if (fMvdBranchName == "Mix")
110  {
111  stripar=(TClonesArray*) ioman->GetObject("MVDHitsStripMix");
112  if(stripar!=0)
113  {
114  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsStripMix"),new GFRecoHitProducer<PndSdsHit,PndSdsRecoHit>(stripar));
115  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "MVDHitsStripMix array found" << std::endl;
116  }
117 
118  pixelar=(TClonesArray*) ioman->GetObject("MVDHitsPixelMix");
119  if(pixelar!=0)
120  {
121  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsPixelMix"),new GFRecoHitProducer<PndSdsHit,PndSdsRecoHit>(pixelar));
122  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "MVDHitsPixelMix array found" << std::endl;
123  }
124  }
125  else
126  {
127  stripar=(TClonesArray*) ioman->GetObject("MVDHitsStrip");
128  if(stripar!=0)
129  {
130  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsStrip"),new GFRecoHitProducer<PndSdsHit,PndSdsRecoHit>(stripar));
131  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "MVDHitsStrip array found" << std::endl;
132  }
133 
134  pixelar=(TClonesArray*) ioman->GetObject("MVDHitsPixel");
135  if(pixelar!=0)
136  {
137  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MVDHitsPixel"),new GFRecoHitProducer<PndSdsHit,PndSdsRecoHit>(pixelar));
138  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "MVDHitsPixel array found" << std::endl;
139  }
140  }
141 
142  TClonesArray* sthit; // TClonesArray *sttr; //[R.K. 01/2017] unused variable?
143 
144  if (fCentralTrackerBranchName == "Mix")
145  {
146  sthit=(TClonesArray*) ioman->GetObject("STTHitMix");
147  if(sthit!=0)
148  {
149  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("STTHitMix"),new PndSttRecoHitProducer<PndSttHit,PndSttRecoHit>(sthit, tubeArray));
150  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "STTHitMix array found" << std::endl;
151  }
152  }
153  else
154  {
155  sthit=(TClonesArray*) ioman->GetObject("STTHit");
156  if(sthit!=0)
157  {
158  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("STTHit"),new PndSttRecoHitProducer<PndSttHit,PndSttRecoHit>(sthit, tubeArray));
159  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "SttHit array found" << std::endl;
160  }
161  }
162 
163  TClonesArray* gemar=(TClonesArray*) ioman->GetObject("GEMHit");
164  if(gemar!=0)
165  {
166  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("GEMHit"),new GFRecoHitProducer<PndGemHit,PndGemRecoHit>(gemar));
167  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "GEMHit array found" << std::endl;
168  }
169 
170  TClonesArray* mdtar=(TClonesArray*) ioman->GetObject("MdtHit");
171  if(mdtar!=0)
172  {
173  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("MdtHit"),new GFRecoHitProducer<PndMdtHit,PndMdtRecoHit>(mdtar));
174  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "MdtHit array found" << std::endl;
175  }
176 
177  TClonesArray* ftsar=(TClonesArray*) ioman->GetObject("FTSHit");
178  if(ftsar!=0)
179  {
180  fTheRecoHitFactory->addProducer(FairRootManager::Instance()->GetBranchId("FTSHit"),new PndFtsRecoHitProducer<PndFtsHit,PndFtsRecoHit>(ftsar, ftsTubeArray));
181  std::cout << "*** PndRecoKalmanFit::Init" << "\t" << "FtsHit array found" << std::endl;
182  }
183 
184  if (fUseGeane)
185  {
186  fPro = new FairGeanePro();
187  if (fVerbose==0) fPro->SetPrintErrors(kFALSE);
188  }
189  else
190  {
191  Error("PndRecoKalmanFit::Init","Only GEANE Propagatio available!!!");
192  return kFALSE;
193  }
194 
196 
198 
199  std::cout << "===PndRecoKalmanFit::Init() finished ===================================================" << std::endl;
200 
201  return kTRUE;
202 }
TString fCentralTrackerBranchName
Name of the TCA for MVD.
Int_t GetGeometryType()
Definition: PndGeoFtsPar.h:37
void addProducer(int detID, GFAbsRecoHitProducer *hitProd)
Register a producer module to the factory.
Short_t fTrackRep
Number of iterations.
FairGeanePro * fPro
Template class for a hit producer module.
Int_t GetGeometryType()
Definition: PndGeoSttPar.h:28
static void quiet(bool b=true)
Definition: GFException.h:93
TString fMvdBranchName
Geane Propagator.
Magnetic field.
Factory object to create RecoHits from digitized and clustered data.
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
Int_t fVerbose
(0) GeaneTrackRep, 1 RKTrackRep
Int_t fNumIt
Flag to use as initial plane the one perpendicular to the track.
Bool_t fUseGeane
Name of the TCA for central tracker.
void init(GFAbsBField *b)
set the magntic field here. Magnetic field classes must be derived from GFAbsBField ...
GFRecoHitFactory * fTheRecoHitFactory
static GFFieldManager * getInstance()
void setNumIterations(Int_t i)
Set number of iterations for Kalman Filter.
Definition: GFKalman.h:86
void PndRecoKalmanFit::SetCentralTrackerBranchName ( const TString name)
inline

Definition at line 48 of file PndRecoKalmanFit.h.

References fCentralTrackerBranchName, and name.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

TString fCentralTrackerBranchName
Name of the TCA for MVD.
TString name
void PndRecoKalmanFit::SetGeane ( Bool_t  opt = kTRUE)
inline

Definition at line 40 of file PndRecoKalmanFit.h.

References fUseGeane.

Referenced by PndRecoMultiKalmanTask::Init(), PndRecoKalmanTask::Init(), and PndPidCorrelator::Init().

40 { fUseGeane = opt; }
Bool_t fUseGeane
Name of the TCA for central tracker.
void PndRecoKalmanFit::SetMvdBranchName ( const TString name)
inline

Definition at line 47 of file PndRecoKalmanFit.h.

References fMvdBranchName, and name.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

47 { fMvdBranchName = name; }
TString fMvdBranchName
Geane Propagator.
TString name
void PndRecoKalmanFit::SetNumIterations ( Int_t  num)
inline

Definition at line 44 of file PndRecoKalmanFit.h.

References fNumIt, and num.

Referenced by PndRecoMultiKalmanTask::Init(), PndRecoKalmanTask::Init(), and PndPidCorrelator::Init().

44 { fNumIt = num; }
int num[96]
Definition: ranlxd.cxx:381
Int_t fNumIt
Flag to use as initial plane the one perpendicular to the track.
void PndRecoKalmanFit::SetPerpPlane ( Bool_t  opt = kTRUE)
inline

Definition at line 43 of file PndRecoKalmanFit.h.

References fPerpPlane.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

43 { fPerpPlane = opt; }
Bool_t fPerpPlane
Distance in [cm] to back-propagate the parameters, negative number means no backpropagation.
void PndRecoKalmanFit::SetPropagateDistance ( Float_t  opt = -1.f)
inline

Definition at line 42 of file PndRecoKalmanFit.h.

References fPropagateDistance.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

42 { fPropagateDistance = opt;}
Float_t fPropagateDistance
Flag to propagate to the interaction point.
void PndRecoKalmanFit::SetPropagateToIP ( Bool_t  opt = kTRUE)
inline

Definition at line 41 of file PndRecoKalmanFit.h.

References fPropagateToIP.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

41 { fPropagateToIP = opt;}
Bool_t fPropagateToIP
Flag to use Geane.
void PndRecoKalmanFit::SetTrackRep ( Int_t  num)
inline

Definition at line 45 of file PndRecoKalmanFit.h.

References fTrackRep, and num.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

45 { fTrackRep = num; }
int num[96]
Definition: ranlxd.cxx:381
Short_t fTrackRep
Number of iterations.
void PndRecoKalmanFit::SetVerbose ( Int_t  verb)
inline

Definition at line 46 of file PndRecoKalmanFit.h.

References fVerbose.

Referenced by PndRecoMultiKalmanTask::Init(), and PndRecoKalmanTask::Init().

46 { fVerbose = verb; }
Int_t fVerbose
(0) GeaneTrackRep, 1 RKTrackRep

Member Data Documentation

TString PndRecoKalmanFit::fCentralTrackerBranchName
private

Name of the TCA for MVD.

Definition at line 64 of file PndRecoKalmanFit.h.

Referenced by Init(), and SetCentralTrackerBranchName().

GFKalman PndRecoKalmanFit::fGenFitter
private

Definition at line 59 of file PndRecoKalmanFit.h.

Referenced by Fit(), and Init().

TString PndRecoKalmanFit::fMvdBranchName
private

Geane Propagator.

Definition at line 63 of file PndRecoKalmanFit.h.

Referenced by Init(), and SetMvdBranchName().

Int_t PndRecoKalmanFit::fNumIt
private

Flag to use as initial plane the one perpendicular to the track.

Definition at line 70 of file PndRecoKalmanFit.h.

Referenced by Init(), and SetNumIterations().

Bool_t PndRecoKalmanFit::fPerpPlane
private

Distance in [cm] to back-propagate the parameters, negative number means no backpropagation.

Definition at line 69 of file PndRecoKalmanFit.h.

Referenced by Fit(), and SetPerpPlane().

FairGeanePro* PndRecoKalmanFit::fPro
private

Definition at line 61 of file PndRecoKalmanFit.h.

Referenced by Fit(), Init(), and ~PndRecoKalmanFit().

Float_t PndRecoKalmanFit::fPropagateDistance
private

Flag to propagate to the interaction point.

Definition at line 68 of file PndRecoKalmanFit.h.

Referenced by Fit(), and SetPropagateDistance().

Bool_t PndRecoKalmanFit::fPropagateToIP
private

Flag to use Geane.

Definition at line 67 of file PndRecoKalmanFit.h.

Referenced by Fit(), and SetPropagateToIP().

GFRecoHitFactory* PndRecoKalmanFit::fTheRecoHitFactory
private

Definition at line 53 of file PndRecoKalmanFit.h.

Referenced by Fit(), GetRecoHitFactory(), and Init().

Short_t PndRecoKalmanFit::fTrackRep
private

Number of iterations.

Definition at line 71 of file PndRecoKalmanFit.h.

Referenced by Fit(), Init(), and SetTrackRep().

Bool_t PndRecoKalmanFit::fUseGeane
private

Name of the TCA for central tracker.

Definition at line 66 of file PndRecoKalmanFit.h.

Referenced by Init(), and SetGeane().

Int_t PndRecoKalmanFit::fVerbose
private

(0) GeaneTrackRep, 1 RKTrackRep

Definition at line 72 of file PndRecoKalmanFit.h.

Referenced by Fit(), Init(), and SetVerbose().


The documentation for this class was generated from the following files: