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

#include <PndHypDKalmanTask.h>

Inheritance diagram for PndHypDKalmanTask:

Public Member Functions

 PndHypDKalmanTask ()
 
 ~PndHypDKalmanTask ()
 
void SetTrackBranchName (const TString &name)
 
void AddHitBranch (unsigned int detId, const TString &m)
 
void SetPersistence (Bool_t opt=kTRUE)
 
void SetSmooth (Bool_t opt=kTRUE)
 
void UseMVD (bool mvd)
 
virtual InitStatus Init ()
 
virtual void Exec (Option_t *opt)
 
void WriteHistograms (const TString &filename)
 

Private Attributes

TString fTrackBranchName
 
std::map< unsigned int, TStringfHitBranchMap
 
TClonesArray * fTrackArray
 
Bool_t fPersistence
 
GFRecoHitFactoryfTheRecoHitFactory
 
FairField * fField
 
TH1D * fPH
 
TH1D * fChi2H
 
TH1D * fXresH
 
TH1D * fYresH
 
TH1D * fXresFitH
 
TH1D * fYresFitH
 
TH1D * fPEnd
 
TH1D * fPull
 
int fTrackcount
 
Bool_t fSmooth
 
int fEvt
 
bool fUseMVD
 

Detailed Description

Definition at line 35 of file PndHypDKalmanTask.h.

Constructor & Destructor Documentation

PndHypDKalmanTask::PndHypDKalmanTask ( )

Definition at line 58 of file PndHypDKalmanTask.cxx.

References fTrackBranchName, and PndGeoHandling::Instance().

59  : FairTask("Kalman Filter"), fPersistence(kFALSE),fSmooth(kFALSE),fUseMVD(false),fEvt(0)
60 {
61  fTrackBranchName = "Track";
63 }
static PndGeoHandling * Instance()
PndHypDKalmanTask::~PndHypDKalmanTask ( )

Definition at line 66 of file PndHypDKalmanTask.cxx.

References fChi2H, and fPH.

67 {
68  if(fPH!=NULL)delete fPH;
69  if(fChi2H!=NULL)delete fChi2H;
70 }

Member Function Documentation

void PndHypDKalmanTask::AddHitBranch ( unsigned int  detId,
const TString m 
)
inline

Definition at line 50 of file PndHypDKalmanTask.h.

References fHitBranchMap, and m.

50 {fHitBranchMap[detId]=m;};
__m128 m
Definition: P4_F32vec4.h:28
std::map< unsigned int, TString > fHitBranchMap
void PndHypDKalmanTask::Exec ( Option_t *  opt)
virtual

Definition at line 159 of file PndHypDKalmanTask.cxx.

References GFTrack::addHitVector(), GFRecoHitFactory::createMany(), fChi2H, fEvt, fPEnd, fPH, fPull, fTheRecoHitFactory, fTrackArray, fTrackcount, GFTrack::getCand(), GFTrack::getCardinalRep(), GFTrack::getChiSqu(), GFAbsTrackRep::getMom(), GFTrack::getMom(), GFTrack::getNumHits(), GFTrack::getNumReps(), GFTrack::getPos(), GFTrackCand::getQoverPseed(), GFTrack::getRedChiSqu(), GFAbsTrackRep::getStatusFlag(), GFTrack::getTrackRep(), p, pos, GFAbsTrackRep::Print(), GFKalman::processTrack(), GFKalman::setNumIterations(), GeaneTrackRep::setPropDir(), trk, and GFException::what().

160 {
161  std::cout<<"PndHypDKalmanTask::Exec Event "<<fEvt++<<std::endl;
162  //Reset output Array
163  //if(_trackArray==0) Fatal("PndHypDKalman::Exec)","No TrackArray");
164  // _trackArray->Delete();
165 
166  Int_t ntracks=fTrackArray->GetEntriesFast();
167 
168 
169 
170 
171  // Fitting ---------------- can go to another task!
172  GFKalman fitter;
173  //fitter.setLazy(1); // tell the fitter to skip hits if error occurs
174  fitter.setNumIterations(1);
175  for(Int_t itr=0;itr<ntracks;++itr){
176  GFTrack* trk=(GFTrack*)fTrackArray->At(itr);
177  GFTrackCand trcnd = trk->getCand();
178 
179  // Load RecoHits
180  try {
182  std::cout<<trk->getNumHits()<<" hits in track "
183  <<itr<<std::endl;
184  }
185 
186  catch(GFException& e) {
187  std::cout << e.what() << std::endl;
188  throw e;
189  }
190 
191  /*
192  // fill tpc residuals
193  std::vector<double> res;
194  trk->getResiduals(2,0,0,res);
195  for(int i=0;i<res.size();++i){
196  _xresH->Fill(res[i]);
197  }
198  res.clear();
199  trk->getResiduals(2,1,0,res);
200  for(int i=0;i<res.size();++i){
201  _yresH->Fill(res[i]);
202  }
203  res.clear();
204  */
205  // Start Fitter
206  try{
207  fitter.processTrack(trk);
208 
209  }
210  catch (GFException e){
211  std::cout<<"*** FITTER EXCEPTION ***"<<std::endl;
212  std::cout<<e.what()<<std::endl;
213 
214  }
215 
216  std::cout<<"SUCESSFULL FIT!"<<std::endl;
217 
218  // Fill some histos after fit
219  if(trk->getTrackRep(0)->getStatusFlag()==0){
220  // propagate backwards
221  GeaneTrackRep* gtrk=dynamic_cast<GeaneTrackRep*>(trk->getTrackRep(0));
222  if(gtrk!=NULL)gtrk->setPropDir(-1);
223 
224  trk->getCardinalRep()->Print();
225 
226  // ------- Propagation to prim vertex ---------
227 
228  GFDetPlane pl(TVector3(0,0,-55.0),TVector3(1,0,0),TVector3(0,1,0));
229 
230  //TVector3 p3=trk->getTrackRep(0)->getMom(pl);
231  //--------------------------------------------
232 
233  double p=trk->getTrackRep(0)->getMom().Mag();
234 
235  std::cout<<" momentum "<<" "<<trk->getMom().Mag()<<std::endl;
236  std::cout<<" "<<" N reps "<<trk->getNumReps()<<" chi2red "<<trk->getRedChiSqu()<<std::endl;
237 
238  //INFO: Changing the EPSIL param to
239  //0.05 (HYPsilicon, HYPdiamond, HYPcarbon--> materials in geo file)
240  //propagation to prim vertex works smoothly for 500 ev.
241  //no dedx modification from media file is needed
242 
243  fPH->Fill(p);
244 
245  fPEnd->Fill(-1/(trcnd.getQoverPseed()));
246 
247  fPull->Fill((-1/(trcnd.getQoverPseed()))-p);
248 
249  TVector3 pos=trk->getPos();
250  //fPEnd->Fill(pos.X(),pos.Z());
251 
252  double chi2=trk->getChiSqu();
253  fChi2H->Fill(chi2);
254  ++fTrackcount;
255 
256 
257 
258  // // fill tpc residuals
259 // trk->getResiduals(2,0,0,res);
260 // for(int i=0;i<res.size();++i){
261 // _xresFitH->Fill(res[i]);
262 // }
263 // res.clear();
264 // trk->getResiduals(2,1,0,res);
265 // for(int i=0;i<res.size();++i){
266 // _yresFitH->Fill(res[i]);
267 // }
268 // res.clear();
269 
270  }
271  }
272 
273  return;
274 }
TVector3 pos
unsigned int getNumHits() const
Definition: GFTrack.h:148
GFRecoHitFactory * fTheRecoHitFactory
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.
double getChiSqu() const
Get chi2.
Definition: GFTrack.h:252
Generic Kalman Filter implementation.
Definition: GFKalman.h:50
TVector3 getPos() const
Get present position.
Definition: GFTrack.h:222
virtual const char * what() const
standard error message handling for exceptions. use like "std::cerr << e.what();" ...
Definition: GFException.cxx:47
void addHitVector(std::vector< GFAbsRecoHit * > hits)
Add collection of hits.
Definition: GFTrack.h:310
TVector3 getMom() const
Get momentum at the present position.
Definition: GFTrack.h:209
Double_t p
Definition: anasim.C:58
double getRedChiSqu() const
Get chi2/NDF.
Definition: GFTrack.h:264
virtual TVector3 getMom(const GFDetPlane &pl)=0
unsigned int getNumReps() const
Get number of track represenatations.
Definition: GFTrack.h:192
bool getStatusFlag()
double getQoverPseed() const
get the seed value for track: qoverp
Definition: GFTrackCand.h:135
Track candidate – a list of cluster indices.
Definition: GFTrackCand.h:55
void setPropDir(int d)
void processTrack(GFTrack *trk)
Performs fit on a GFTrack.
Definition: GFKalman.cxx:40
GFAbsTrackRep * getTrackRep(int id) const
Accessor for track representations.
Definition: GFTrack.h:186
GFTrack * trk
Definition: checkgenfit.C:13
GFAbsTrackRep * getCardinalRep() const
Get cardinal track representation.
Definition: GFTrack.h:202
virtual void Print() const
const GFTrackCand & getCand() const
Definition: GFTrack.h:142
Exception class for error handling in GENFIT (provides storage for diagnostic information) ...
Definition: GFException.h:45
TClonesArray * fTrackArray
void setNumIterations(Int_t i)
Set number of iterations for Kalman Filter.
Definition: GFKalman.h:86
InitStatus PndHypDKalmanTask::Init ( )
virtual

Definition at line 73 of file PndHypDKalmanTask.cxx.

References GFRecoHitFactory::addProducer(), fChi2H, fHitBranchMap, fPEnd, fPH, fPull, fTheRecoHitFactory, fTrackArray, fTrackBranchName, fUseMVD, fXresFitH, fXresH, fYresFitH, fYresH, and gGeoManager.

74 {
75  //Get ROOT Manager
76  FairRootManager* ioman= FairRootManager::Instance();
77 
78  if(ioman==0)
79  {
80  Error("PndHypDKalmanTask::Init","RootManager not instantiated!");
81  return kERROR;
82  }
83 
84  // Get input collection
85  fTrackArray=(TClonesArray*) ioman->GetObject(fTrackBranchName);
86 
87  //if(_trackArray==0)
88  //{
89  // Error("PndHypDKalmanTask::Init","track-array not found!");
90  //return kERROR;
91  //}
92 
93  // create and register output array
94  //_trackArray = new TClonesArray("GFTrack");
95  //ioman->Register("TrackPreFit","GenFit",_trackArray,_persistence);
96 
97  // Build hit factory -----------------------------
99 
100  std::map<unsigned int,TString>::iterator iter=fHitBranchMap.begin();
101  while(iter!=fHitBranchMap.end()){
102 
103  TClonesArray* har=(TClonesArray*) ioman->GetObject("HypHit");
104 
105  //cout<<" "<<iter->second<<endl;
106 
107  if(har==0){
108 
109  Error("PndHypKalmanTask::Init","Hit array not found");
110  }
111  else{
112 
113  if(iter->first==2)fTheRecoHitFactory->addProducer(iter->first,new GFRecoHitProducer<PndHypHit,PndHypRecoSPHit>(har));
114 
115  }
116 
117  if(fUseMVD==true){
118  TClonesArray* sar=(TClonesArray*) ioman->GetObject("MVDHit");
119  if(sar==0){
120 
121  Error("PndHypKalmanTask::Init","Hit array not found");
122  }
123  else{
124 
125  if(iter->first==3)fTheRecoHitFactory->addProducer(iter->first,new GFRecoHitProducer<PndSdsHit,PndSdsRecoHit>(sar));
126  }
127  }
128 
129  // TClonesArray* ar=(TClonesArray*) ioman->GetObject(iter->second);
130  // if(ar==0){
131  // Error("PndHypDKalmanTask::Init","point-array %s not found!",iter->second.Data());
132  // }
133  // else{
134  // // the next lines is not general because it will work only for CmMCPoints!
135  // fTheRecoHitFactory->addProducer(iter->first,new GFRecoHitProducer<PndHypHit,PndHypRecoSPHit>(ar));
136  // //FairHit,PndHypDSPHit>(ar));PndHypPoint,PndHypRecoHit
137  // }
138 
139 
140  ++iter;
141  }//end loops over hit types
142 
143  // setup histograms
144  fPH=new TH1D("pH","p",500,0.02,0.7);
145  fChi2H=new TH1D("chi2H","chi2",100,0,20);
146  fXresH=new TH1D("xres","xres",100,-5,5);
147  fYresH=new TH1D("yres","yres",100,-5,5);
148  fXresFitH=new TH1D("xresfit","xres after fit",100,-5,5);
149  fYresFitH=new TH1D("yresfit","yres after fit",100,-5,5);
150  fPEnd=new TH1D("pPre","Endpoint",500,0.02,0.7);
151  fPull=new TH1D("pPull","Pull",500,-0.3,0.3);
152  std::cout << "-I- gGeoManager = "<<gGeoManager << std::endl;
153 
154  return kSUCCESS;
155 }
GFRecoHitFactory * fTheRecoHitFactory
void addProducer(int detID, GFAbsRecoHitProducer *hitProd)
Register a producer module to the factory.
Template class for a hit producer module.
TGeoManager * gGeoManager
std::map< unsigned int, TString > fHitBranchMap
Factory object to create RecoHits from digitized and clustered data.
TClonesArray * fTrackArray
void PndHypDKalmanTask::SetPersistence ( Bool_t  opt = kTRUE)
inline

Definition at line 51 of file PndHypDKalmanTask.h.

References fPersistence.

51 {fPersistence=opt;}
void PndHypDKalmanTask::SetSmooth ( Bool_t  opt = kTRUE)
inline

Definition at line 52 of file PndHypDKalmanTask.h.

References fSmooth.

52 {fSmooth=opt;}
void PndHypDKalmanTask::SetTrackBranchName ( const TString name)
inline

Definition at line 49 of file PndHypDKalmanTask.h.

References fTrackBranchName, and name.

TString name
void PndHypDKalmanTask::UseMVD ( bool  mvd)
inline

Definition at line 53 of file PndHypDKalmanTask.h.

References fUseMVD.

53 {fUseMVD=mvd;}
void PndHypDKalmanTask::WriteHistograms ( const TString filename)

Definition at line 277 of file PndHypDKalmanTask.cxx.

References fChi2H, file, fPEnd, fPH, fPull, fXresFitH, fXresH, fYresFitH, and fYresH.

277  {
278  TFile* file = FairRootManager::Instance()->GetOutFile();
279  file->cd();
280  file->mkdir("PndHypDKalman");
281  file->cd("PndHypDKalman");
282 
283  // TFile* file = new TFile(filename,"UPDATE");
284  // if(file->cd("Kalman")) file->Delete("Kalman;*");
285  // file->mkdir("Kalman");
286  // file->cd("Kalman");
287 
288  fPH->Write();
289  delete fPH;
290  fPH=NULL;
291 
292  fChi2H->Write();
293  delete fChi2H;
294  fChi2H=NULL;
295 
296  fXresH->Write();
297  delete fXresH;
298  fXresH=NULL;
299 
300  fYresH->Write();
301  delete fYresH;
302  fYresH=NULL;
303 
304  fXresFitH->Write();
305  delete fXresFitH;
306  fXresFitH=NULL;
307 
308  fYresFitH->Write();
309  delete fYresFitH;
310  fYresFitH=NULL;
311 
312  fPEnd->Write();
313  delete fPEnd;
314  fPEnd=NULL;
315 
316  fPull->Write();
317  delete fPull;
318  fPull=NULL;
319 
320  file->Close();delete file;
321 
322 }
TFile * file

Member Data Documentation

TH1D* PndHypDKalmanTask::fChi2H
private

Definition at line 76 of file PndHypDKalmanTask.h.

Referenced by Exec(), Init(), WriteHistograms(), and ~PndHypDKalmanTask().

int PndHypDKalmanTask::fEvt
private

Definition at line 86 of file PndHypDKalmanTask.h.

Referenced by Exec().

FairField* PndHypDKalmanTask::fField
private

Definition at line 73 of file PndHypDKalmanTask.h.

std::map<unsigned int,TString> PndHypDKalmanTask::fHitBranchMap
private

Definition at line 67 of file PndHypDKalmanTask.h.

Referenced by AddHitBranch(), and Init().

TH1D* PndHypDKalmanTask::fPEnd
private

Definition at line 81 of file PndHypDKalmanTask.h.

Referenced by Exec(), Init(), and WriteHistograms().

Bool_t PndHypDKalmanTask::fPersistence
private

Definition at line 70 of file PndHypDKalmanTask.h.

Referenced by SetPersistence().

TH1D* PndHypDKalmanTask::fPH
private

Definition at line 75 of file PndHypDKalmanTask.h.

Referenced by Exec(), Init(), WriteHistograms(), and ~PndHypDKalmanTask().

TH1D* PndHypDKalmanTask::fPull
private

Definition at line 82 of file PndHypDKalmanTask.h.

Referenced by Exec(), Init(), and WriteHistograms().

Bool_t PndHypDKalmanTask::fSmooth
private

Definition at line 85 of file PndHypDKalmanTask.h.

Referenced by SetSmooth().

GFRecoHitFactory* PndHypDKalmanTask::fTheRecoHitFactory
private

Definition at line 72 of file PndHypDKalmanTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndHypDKalmanTask::fTrackArray
private

Definition at line 68 of file PndHypDKalmanTask.h.

Referenced by Exec(), and Init().

TString PndHypDKalmanTask::fTrackBranchName
private

Definition at line 66 of file PndHypDKalmanTask.h.

Referenced by Init(), PndHypDKalmanTask(), and SetTrackBranchName().

int PndHypDKalmanTask::fTrackcount
private

Definition at line 84 of file PndHypDKalmanTask.h.

Referenced by Exec().

bool PndHypDKalmanTask::fUseMVD
private

Definition at line 87 of file PndHypDKalmanTask.h.

Referenced by Init(), and UseMVD().

TH1D* PndHypDKalmanTask::fXresFitH
private

Definition at line 79 of file PndHypDKalmanTask.h.

Referenced by Init(), and WriteHistograms().

TH1D* PndHypDKalmanTask::fXresH
private

Definition at line 77 of file PndHypDKalmanTask.h.

Referenced by Init(), and WriteHistograms().

TH1D* PndHypDKalmanTask::fYresFitH
private

Definition at line 80 of file PndHypDKalmanTask.h.

Referenced by Init(), and WriteHistograms().

TH1D* PndHypDKalmanTask::fYresH
private

Definition at line 78 of file PndHypDKalmanTask.h.

Referenced by Init(), and WriteHistograms().


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