FairRoot/PandaRoot
Functions
rho/tut_ana_mclist.C File Reference

Go to the source code of this file.

Functions

void tut_ana_mclist (int nevts=0, TString prefix="signal")
 

Function Documentation

void tut_ana_mclist ( int  nevts = 0,
TString  prefix = "signal" 
)

Definition at line 8 of file rho/tut_ana_mclist.C.

References PndAnalysis::FillList(), fRun, PndAnalysis::GetEntries(), PndAnalysis::GetEvent(), RhoCandList::GetLength(), RhoCandidate::GetTrackNumber(), i, PndAnalysis::PndAnalysis(), rtdb, and TString.

9 {
10  // *** some variables
11  int i=0,j=0, k=0, l=0;
12  gStyle->SetOptFit(1011);
13 
14  // *** the output file for FairRunAna
15  TString OutFile="out_dummy.root";
16 
17  // *** the files coming from the simulation
18  TString inPidFile = prefix+"_pid.root"; // this file contains the PndPidCandidates and McTruth
19  TString inParFile = prefix+"_par.root";
20 
21  // *** PID table with selection thresholds; can be modified by the user
22  TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
23 
24  // *** initialization
25  FairLogger::GetLogger()->SetLogToFile(kFALSE);
26  FairRunAna* fRun = new FairRunAna();
27  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
28  fRun->SetSource(new FairFileSource(inPidFile));
29 
30  // *** setup parameter database
31  FairParRootFileIo* parIO = new FairParRootFileIo();
32  parIO->open(inParFile);
33  FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
34  parIOPid->open(pidParFile.Data(),"in");
35 
36  rtdb->setFirstInput(parIO);
37  rtdb->setSecondInput(parIOPid);
38  rtdb->setOutput(parIO);
39 
40  fRun->SetOutputFile(OutFile);
41  fRun->Init();
42 
43  //
44  // Now the analysis stuff comes...
45  //
46 
47 
48  // *** the data reader object
49  PndAnalysis* theAnalysis = new PndAnalysis();
50  if (nevts==0) nevts= theAnalysis->GetEntries();
51 
52  // *** RhoCandLists for the analysis
53  RhoCandList mctruth;
54 
55  // ***
56  // the event loop
57  // ***
58  while (theAnalysis->GetEvent() && i++<nevts)
59  {
60  cout<<"****** Evt " << i << endl;
61 
62  // *** the MC Truth objects
63  theAnalysis->FillList(mctruth,"McTruth");
64 
65  //
66  // Print MC Truth list with mother-daughter relations
67  //
68  for (j=0;j<mctruth.GetLength();++j)
69  {
70  RhoCandidate *mcmother = mctruth[j]->TheMother();
71 
72  int muid = -1;
73  if (mcmother) muid = mcmother->GetTrackNumber();
74 
75  cout << "Track "<< mctruth[j]->GetTrackNumber()<<" (PDG:"<<mctruth[j]->PdgCode() <<") has mother "<<muid;
76 
77  if (mctruth[j]->NDaughters()>0) cout <<" and daughter(s) ";
78  for (k=0;k<mctruth[j]->NDaughters();++k) cout <<mctruth[j]->Daughter(k)->GetTrackNumber()<<" ";
79 
80  cout<<endl;
81  }
82  cout <<endl;
83  }
84 
85 }
Int_t GetEntries()
Int_t GetTrackNumber() const
Definition: RhoCandidate.h:417
Int_t i
Definition: run_full.C:25
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
FairRunAna * fRun
Definition: hit_dirc.C:58
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
Int_t GetEvent(Int_t n=-1)