FairRoot/PandaRoot
anaDMesonsCharged.C
Go to the documentation of this file.
1 int anaDMesonsCharged(TString filename="data/dpdm", int nevts=0)
2 {
3  TStopwatch timer;
4  timer.Start();
5  int verbose=0;
6  int canvasplot=1; //whether to plot into a ps, too.
7  int canvasdraw=1;
8  gSystem->Load("libEGPythia6"); // needed for TDatabasePDG
9  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
10  gSystem->Load("libRho");
11  gSystem->Load("libPid");
12  gSystem->Load("libAnalysisTools");
13 
14  gROOT->LoadMacro("plothistos.C");
15  //gROOT->LoadMacro("$VMCWORKDIR/macro/mvd/Tools.C");
16  //LoadPandaStyle();
17 
18  int bins=100;
19  double min=0,max=2.5; // GeV/c
20  double min2=0.13,max2=0.15;
21  double min3=1.,max3=1.5;
22  double max4=0.3;//cm
23  double maxchi = 50;
24 
25  TList* histolist = new TList();
26 
27  //TH1F *hChisquare=new TH1F("chisq","#Chi^{2}",bins,0,maxchi);
30  //TH2F *hChisquareMom=new TH2F("chisqmom",";p / GeV/c;#Chi^{2}",bins,min,max,bins,0,maxchi);
31  //TH2F *hChisquareTheta=new TH2F("chisqthe",";#theta / rad;#Chi^{2}",bins,0,3.5,bins,0,maxchi);
32  //TH2F *hMomMass=new TH2F("mommass","Mom vs mass;M*Q / e*GeV/c^{2};p / GeV/c",bins,-max2,max2,bins,min,max);
33 
34  TH1F *hMompip=new TH1F("momppip","P(#pi+);p / GeV/c",bins,min,max);
35  TH1F *hMompim=new TH1F("momppim","P(#pi-);p / GeV/c",bins,min,max);
36  TH1F *hMomKp=new TH1F("mompkp","P(K+);p / GeV/c",bins,min,max);
37  TH1F *hMomKm=new TH1F("mompkm","P(K-);p / GeV/c",bins,min,max);
38  TH1I *hTFStat=new TH1I("tfstat","Track Fit Status Flag",10,-5,5);
39 
40  TH1F *hDpMass=new TH1F("dpmass","Mass D+;m / (GeV/c^{2})",bins,0,5);
41  TH1F *hDmMass=new TH1F("dmmass","Mass D-;m / (GeV/c^{2})",bins,0,5);
42  TH1F *hBeamMass=new TH1F("beammass","Mass Beam;m / (GeV/c^{2})",bins,0,5);
43 
44  TH1F *hVtxX=new TH1F("hvtxx","Vertex after fit (D+ and D-);x/cm",bins,-0.2,0.2);
45  TH1F *hVtxY=new TH1F("hvtxy","Vertex after fit (D+ and D-);y/cm",bins,-0.2,0.2);
46  TH1F *hVtxZ=new TH1F("hvtxz","Vertex after fit (D+ and D-);z/cm",bins,-0.2,0.2);
47  TH1F *hVtxR=new TH1F("hvtxr","Vertex after fit (D+ and D-);|r|/cm",bins,0,20);
48  TH1F *hVtxChi2=new TH1F("hvtxchi2","Vertex Fit #Chi^{2};#Chi^{2}",bins,0,10);
49 
50  // Add histogram pointers to a list for iterative saving & plotting
51  //---newpage
52  histolist->Add(hMompip);
53  histolist->Add(hMompim);
54  histolist->Add(hMomKp);
55  histolist->Add(hMomKm);
56  histolist->Add(hDpMass);
57  histolist->Add(hDmMass);
58  histolist->Add(hBeamMass);
59  histolist->Add(hTFStat);
60  histolist->Add(hVtxX);
61  histolist->Add(hVtxY);
62  histolist->Add(hVtxZ);
63  histolist->Add(hVtxR);
64  histolist->Add(hVtxChi2);
65 
66  // Set up the Analysis
67  TString histoFile = InitRun(filename);
68  PndAnalysis* theAnalysis = new PndAnalysis();
69 
70  // the candidates lists we need
71  TCandList pionp, pionm, kaonp, kaonm, dplus, dminus ,mctracks, beam;
72  TCandidate *tmpcand=0, *recocand=0, *mccand=0;
73 
74  int evts = theAnalysis->GetEntries();
75  if (nevts>0 && nevts<evts) evts=nevts;
76 
77  double tfcut=-1; // trackfit status cut
78 
79  TLorentzVector tmpMom;
80  TVector3 p1,p2,backboost,vertex;
81  double mass=0,momentum=0,theta=0,chi2,energy=0;
82  int ievt=0,j=0,mcidx=-1,fitstat=0;
83  bool tester=false;
84 
85  TDatabasePDG *pdgDatabase = TRho::Instance()->GetPDG(); // Access particle DB
86  //TParticlePDG* pTypePi0=pdgDatabase->GetParticle(111); // Pi0 PDG entry
87  // some mass selectors
88  //TPidMassSelector *piMassSel=new TPidMassSelector("pi0",pTypePi0->Mass(),0.03);
89 
90  // *************
91  // this is the loop through the events ... as simple as this...
92  // ****************
93  std::cout<<" ***** Start the eventloop with "<<evts<<" events. *****"<<std::endl;
94  while (theAnalysis->GetEvent() && ievt++<evts)
95  {
96  theAnalysis->FillList(mctracks,"McTruth");
97  theAnalysis->FillList(pionp,"PionLoosePlus");
98  theAnalysis->FillList(pionm,"PionLooseMinus");
99  theAnalysis->FillList(kaonp,"KaonLoosePlus");
100  theAnalysis->FillList(kaonm,"KaonLooseMinus");
101  if(verbose>0) printf("Event %i: TCandList sizes: pionp:%i, pionm:%i, kaonp:%i, kaonm:%i, mctracks:%i\n",ievt,pionp.GetLength(),pionm.GetLength(),kaonp.GetLength(),kaonm.GetLength(),mctracks.GetLength());
102 
103  for (j=0;j<pionp.GetLength();++j) {
104  recocand=pionp.Get(j);
105  if(!recocand) {if(verbose>0) cout<<"############## pionp\n#############"<<endl;continue;}
106  if(verbose>1) std::cout<<" Pion+: "<<*recocand<<std::endl;
107  fitstat=recocand->GetMicroCandidate().GetFitStatus();
108  hTFStat->Fill(fitstat);
109  if(fitstat<tfcut) {
110  pionm.Remove(*recocand);
111  continue;
112  }
113  tmpMom=recocand->P4();
114  momentum=tmpMom.P();
115  hMompip->Fill(momentum);
116  }
117  for (j=0;j<pionm.GetLength();++j) {
118  recocand=pionm.Get(j);
119  if(!recocand) {if(verbose>0) cout<<"############## pionm\n#############"<<endl;continue;}
120  if(verbose>1) std::cout<<" Pion-: "<<*recocand<<std::endl;
121  fitstat=recocand->GetMicroCandidate().GetFitStatus();
122  hTFStat->Fill(fitstat);
123  if(fitstat<tfcut) {
124  pionm.Remove(*recocand);
125  continue;
126  }
127  tmpMom=recocand->P4();
128  momentum=tmpMom.P();
129  hMompim->Fill(momentum);
130  }
131  for (j=0;j<kaonp.GetLength();++j) {
132  recocand=kaonp.Get(j);
133  if(!recocand) {if(verbose>0) cout<<"############## kaonp\n#############"<<endl;continue;}
134  if(verbose>1) std::cout<<" K+: "<<*recocand<<std::endl;
135  fitstat=recocand->GetMicroCandidate().GetFitStatus();
136  hTFStat->Fill(fitstat);
137  if(fitstat<tfcut) {
138  kaonp.Remove(*recocand);
139  continue;
140  }
141  tmpMom=recocand->P4();
142  momentum=tmpMom.P();
143  hMomKp->Fill(momentum);
144  }
145  for (j=0;j<kaonm.GetLength();++j) {
146  recocand=kaonm.Get(j);
147  if(!recocand) {if(verbose>0) cout<<"############## kaonm\n#############"<<endl;continue;}
148  if(verbose>1) std::cout<<" K-: "<<*recocand<<std::endl;
149  fitstat=recocand->GetMicroCandidate().GetFitStatus();
150  hTFStat->Fill(fitstat);
151  if(fitstat<tfcut) {
152  kaonm.Remove(*recocand);
153  continue;
154  }
155  tmpMom=recocand->P4();
156  momentum=tmpMom.P();
157  hMomKm->Fill(momentum);
158  }
159 
160  // for (j=0;j<pionp.GetLength();++j) {
161  // recocand = pionp.Get(j);
162  // if(!recocand) {if(verbose>0) cout<<"############## pi+\n#############"<<endl;continue;}
163  // chi2=recocand->GetMicroCandidate().GetChiSquared();
164  // if(0==chi2) continue;
165  // if(verbose>1) std::cout<<" Pion+: "<<*recocand<<std::endl;
166  // if(verbose>1) std::cout<<" recocand->Uid()="<<recocand->Uid()<<std::endl;
167  // tmpMom=recocand->P4();
168  // mass=tmpMom.M(); momentum=tmpMom.P(); theta=tmpMom.Theta();
169  // mcidx=recocand->GetMicroCandidate().GetMcIndex();
170  // mccand=mctracks.Get(mcidx);
171  // if(verbose>2) printf("Trying mccand no %i at position %p from pi+ candidate %i at %p.\n",mcidx,mccand,j,recocand);
172  // if(!mccand || mcidx<0) {if(verbose>0) printf("did not find mccand for pi+, mcid=%i\n",mcidx); continue;}
173  // }
174 
175  dplus.Combine(kaonm,pionp,pionp);
176  //Vertex Fit. This changes daughters, too!
177  for (j=0;j<dplus.GetLength();++j) {
178  recocand=dplus.Get(j);
179  if(!recocand) {if(verbose>0) cout<<"############## Combine dplus\n#############"<<endl;continue;}
180  if(verbose>1) std::cout<<" recocand->Uid()="<<recocand->Uid()<<std::endl;
181  if(verbose>1) recocand->Print();
182  //PndVtxFitter fitter(*recocand); //D. Mishra 2008
183  //PndChiVtxFitter fitter(*recocand); //V. Jha 2010 defunc-gives zeros?
184  RhoKinVtxFitter fitter(*recocand); //V. Jha 2010 to be tested
185  //RhoKinFitter fitter(*recocand); //V. Jha 2010
186  //fitter.SetVerbose();
187  fitter.Fit();
188  recocand=fitter.FittedCand(*recocand);
189  if(!recocand) {if(verbose>0) cout<<"############## Fit dplus\n#############"<<endl;continue;}
190  if(verbose>1) std::cout<<" recocand->Uid()="<<recocand->Uid()<<std::endl;
191  if(verbose>1) recocand->Print();
192  // plot
193  tmpMom=recocand->P4();
194  mass=tmpMom.M();
195  hDpMass->Fill(mass);
196  vertex=recocand->GetPosition();
197  if(verbose>2) std::cout<<"Vertex of the fit cand is:"<<std::endl;
198  if(verbose>2) vertex.Print();
199  hVtxX->Fill(vertex.X());
200  hVtxY->Fill(vertex.Y());
201  hVtxZ->Fill(vertex.Z());
202  hVtxR->Fill(vertex.Mag());
203  hVtxChi2->Fill(fitter.GlobalChi2());
204  }
205  dminus.Combine(kaonp,pionm,pionm);
206  //Vertex Fit. This changes daughters, too!
207  for (j=0;j<dminus.GetLength();++j) {
208  recocand=dminus.Get(j);
209  if(!recocand) {if(verbose>0) cout<<"############## Combine dminus\n#############"<<endl;continue;}
210  if(verbose>1) std::cout<<" recocand->Uid()="<<recocand->Uid()<<std::endl;
211  if(verbose>1) recocand->Print();
212  //PndVtxFitter fitter(*recocand); //D. Mishra 2008
213  //PndChiVtxFitter fitter(*recocand); //V. Jha 2010 defunc-gives zeros?
214  RhoKinVtxFitter fitter(*recocand); //V. Jha 2010 to be tested
215  //RhoKinFitter fitter(*recocand); //V. Jha 2010
216  //fitter.SetVerbose();
217  fitter.Fit();
218  recocand=fitter.FittedCand(*recocand);
219  if(!recocand) {if(verbose>0) cout<<"############## Fit dminus\n#############"<<endl;continue;}
220  if(verbose>1) std::cout<<" recocand->Uid()="<<recocand->Uid()<<std::endl;
221  if(verbose>1) recocand->Print();
222  // plot
223  tmpMom=recocand->P4();
224  mass=tmpMom.M();
225  hDmMass->Fill(mass);
226  vertex=recocand->GetPosition();
227  if(verbose>2) std::cout<<"Vertex of the fit cand is:"<<std::endl;
228  if(verbose>2) vertex.Print();
229  hVtxX->Fill(vertex.X());
230  hVtxY->Fill(vertex.Y());
231  hVtxZ->Fill(vertex.Z());
232  hVtxR->Fill(vertex.Mag());
233  hVtxChi2->Fill(fitter.GlobalChi2());
234  }
235 
236  beam.Combine(dplus,dminus);
237  for (j=0;j<beam.GetLength();++j) {
238  recocand=beam.Get(j);
239  if(!recocand) {if(verbose>0) cout<<"############## beam\n#############"<<endl;continue;}
240  if(verbose>1) std::cout<<" beam: "<<*recocand<<std::endl;
241  tmpMom=recocand->P4();
242  mass=tmpMom.M(); momentum=tmpMom.P(); theta=tmpMom.Theta();
243  hBeamMass->Fill(mass);
244  }
245  if(verbose>0) cout << " ----- event done -----"<<endl;
246  } // end eventloop
247 
248  // Save all histograms to a file
249  TFile* histofile = new TFile(histoFile.Data(),"RECREATE");
250  histofile->cd();
251 
252  TListIter iter(histolist);
253  TH1* tmph=NULL;
254  while( tmph=(TH1*)iter() ) tmph->Write();
255  histofile->Close();
256 
257  if(canvasplot){
258  plothistos(histoFile);
259  }
260 
261  if(canvasdraw){
262  TCanvas* can = new TCanvas();
263  //can->Divide(2,2);
264  //can->cd(1);
265  gPad->SetLogy();
266  hVtxR->Draw();
267  hVtxR->Fit("exp");
268  //can->cd(2);
269  //can->cd(3);
270  //can->cd(4);
271  }
272 
273  // now take the time...
274  timer.Stop();
275  Double_t rtime = timer.RealTime();
276  Double_t ctime = timer.CpuTime();
277  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
278  return 0;
279 }
280 
281 
283 {
284  FairRunAna* fRun = new FairRunAna();
285  FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
286 
287  filename.ReplaceAll("_sim.root",".root");
288  PndFileNameCreator namecreator(filename.Data());
292  TString tracksFile = namecreator.GetCustomFileName("pidtracks");
293  TString histoFile = namecreator.GetCustomFileName("histos");
294  TString evrdummy = namecreator.GetCustomFileName("evrdummy");
295  cout<<"simFile="<<simFile.Data()<<endl;
296  cout<<"parFile="<<parFile.Data()<<endl;
297  cout<<"recoFile="<<recoFile.Data()<<endl;
298  cout<<"trkFile="<<tracksFile.Data()<<endl;
299  fRun->SetInputFile(simFile);
300  fRun->AddFriend(recoFile);//,"rec");
301  fRun->AddFriend(tracksFile);//,"trk");
302  FairParRootFileIo* parIO = new FairParRootFileIo();
303  parIO->open(parFile.Data());
304  rtdb->setFirstInput(parIO);
305  rtdb->setOutput(parIO);
306 
307  TString outFile = evrdummy;
308  TString sysFile = gSystem->Getenv("VMCWORKDIR");
309 
310  fRun->SetOutputFile(outFile.Data());
311  FairGeane* geane = new FairGeane();
312  fRun->AddTask(geane);
313  fRun->Init();
314  return histoFile;
315 }
316 
317 
318 
Int_t GetEntries()
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
std::string GetParFileName(std::string addon="", bool cut=false)
int anaDMesonsCharged(TString filename="data/dpdm", int nevts=0)
TString outFile
Definition: hit_dirc.C:17
TString InitRun(TString filename)
#define verbose
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
std::string GetSimFileName(std::string addon="", bool cut=false)
FairRunAna * fRun
Definition: hit_dirc.C:58
TString simFile
Definition: bump_emc.C:11
std::string GetRecoFileName(std::string addon="", bool cut=false)
TString sysFile
A simple class which adds the corresponding file extensions to a given base class.
Double_t
TString parFile
Definition: hit_dirc.C:14
PndAnalysis(TString tname1="", TString tname2="", TString algnamec="PidAlgoIdealCharged", TString algnamen="PidAlgoIdealNeutral")
Definition: PndAnalysis.cxx:48
TStopwatch timer
Definition: hit_dirc.C:51
std::string GetCustomFileName(std::string ext, std::string addon="", bool cut=false)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
PndFileNameCreator namecreator("../data/Lars/MvdDtsSim.root")
Double_t ctime
Definition: hit_dirc.C:114
TPad * p2
Definition: hist-t7.C:117
std::string recoFile
TPad * p1
Definition: hist-t7.C:116
Int_t GetEvent(Int_t n=-1)
Double_t rtime
Definition: hit_dirc.C:113
plothistos(TString filename="histos.root")
Definition: plothistos.C:1
Double_t energy
Definition: plot_dirc.C:15
const string filename