FairRoot/PandaRoot
Functions
rho/tut_ana_fit.C File Reference

Go to the source code of this file.

Functions

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

Function Documentation

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

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

References RhoKinFitter::Add4MomConstraint(), RhoKinFitter::AddMassConstraint(), RhoCandList::Combine(), PndAnalysis::FillList(), RhoKinFitter::Fit(), RhoFitterBase::Fit(), fRun, RhoFitterBase::GetChi2(), PndAnalysis::GetEntries(), PndAnalysis::GetEvent(), RhoCandList::GetLength(), RhoFitterBase::GetProb(), hvpos, i, RhoCandidate::M(), muminus, muplus, out, piminus, piplus, PndAnalysis::PndAnalysis(), RhoCandidate::Pos(), rtdb, RhoCandList::Select(), 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  // *** create an output file for all histograms
44  TFile *out = TFile::Open(prefix+"_ana_fit.root","RECREATE");
45 
46  // *** create some histograms
47  TH1F *hjpsim_all = new TH1F("hjpsim_all","J/#psi mass",200,0,4.5);
48  TH1F *hpsim_all = new TH1F("hpsim_all","#psi(2S) mass",200,0,5);
49 
50  TH1F *hjpsim_vf = new TH1F("hjpsim_vf", "J/#psi mass (vertex fit)",200,0,4.5);
51  TH1F *hjpsim_4cf = new TH1F("hjpsim_4cf","J/#psi mass (4C fit)",200,0,4.5);
52  TH1F *hjpsim_mcf = new TH1F("hjpsim_mcf","J/#psi mass (mass constraint fit)",200,0,4.5);
53 
54  TH1F *hjpsi_chi2_vf = new TH1F("hjpsi_chi2_vf", "J/#psi: #chi^{2} vertex fit",100,0,10);
55  TH1F *hpsi_chi2_4c = new TH1F("hpsi_chi2_4c", "#psi(2S): #chi^{2} 4C fit",100,0,250);
56  TH1F *hjpsi_chi2_mf = new TH1F("hjpsi_chi2_mf", "J/#psi: #chi^{2} mass fit",100,0,10);
57 
58  TH1F *hjpsi_prob_vf = new TH1F("hjpsi_prob_vf", "J/#psi: Prob vertex fit",100,0,1);
59  TH1F *hpsi_prob_4c = new TH1F("hpsi_prob_4c", "#psi(2S): Prob 4C fit",100,0,1);
60  TH1F *hjpsi_prob_mf = new TH1F("hjpsi_prob_mf", "J/#psi: Prob mass fit",100,0,1);
61 
62  TH2F *hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-2,2,100,-2,2);
63 
64  //
65  // Now the analysis stuff comes...
66  //
67 
68 
69  // *** the data reader object
70  PndAnalysis* theAnalysis = new PndAnalysis();
71  if (nevts==0) nevts= theAnalysis->GetEntries();
72 
73  // *** RhoCandLists for the analysis
74  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
75 
76  // *** Mass selector for the jpsi cands
77  double m0_jpsi = TDatabasePDG::Instance()->GetParticle("J/psi")->Mass(); // Get nominal PDG mass of the J/psi
78  RhoMassParticleSelector *jpsiMassSel=new RhoMassParticleSelector("jpsi",m0_jpsi,1.0);
79 
80  // *** the lorentz vector of the initial psi(2S), needed by 4C fitter
81  TLorentzVector ini(0, 0, 6.231552, 7.240065);
82 
83  // ***
84  // the event loop
85  // ***
86  while (theAnalysis->GetEvent() && i++<nevts)
87  {
88  if ((i%100)==0) cout<<"evt " << i << endl;
89 
90  // *** Select with no PID info ('All'); type and mass are set
91  theAnalysis->FillList(muplus, "MuonAllPlus");
92  theAnalysis->FillList(muminus, "MuonAllMinus");
93  theAnalysis->FillList(piplus, "PionAllPlus");
94  theAnalysis->FillList(piminus, "PionAllMinus");
95 
96  // *** combinatorics for J/psi -> mu+ mu-
97  jpsi.Combine(muplus, muminus);
98 
99  // ***
100  // *** do VERTEX FIT (J/psi)
101  // ***
102  for (j=0;j<jpsi.GetLength();++j)
103  {
104  hjpsim_all->Fill(jpsi[j]->M()); // fill histo for all J/psi candidates
105 
106  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
107  vtxfitter.Fit();
108 
109  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
110  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
111  hjpsi_chi2_vf->Fill(chi2_vtx);
112  hjpsi_prob_vf->Fill(prob_vtx);
113 
114  if ( prob_vtx > 0.01 ) // when good enough, fill some histos
115  {
116  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
117  TVector3 jVtx=jfit->Pos(); // and the decay vertex position
118 
119  hjpsim_vf->Fill(jfit->M());
120  hvpos->Fill(jVtx.X(),jVtx.Y());
121  }
122  }
123 
124  // *** some rough mass selection
125  jpsi.Select(jpsiMassSel);
126 
127  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
128  psi2s.Combine(jpsi, piplus, piminus);
129 
130  // ***
131  // *** do 4C FIT (initial psi(2S) system)
132  // ***
133  for (j=0;j<psi2s.GetLength();++j)
134  {
135  hpsim_all->Fill(psi2s[j]->M()); // fill histo for all psi(2S) candidates
136 
137  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
138  fitter.Add4MomConstraint(ini); // set 4 constraint
139  fitter.Fit(); // do fit
140 
141  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
142  double prob_4c = fitter.GetProb(); // access probability of fit
143  hpsi_chi2_4c->Fill(chi2_4c);
144  hpsi_prob_4c->Fill(prob_4c);
145 
146  if ( prob_4c > 0.01 ) // when good enough, fill some histo
147  {
148  RhoCandidate *jfit = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
149 
150  hjpsim_4cf->Fill(jfit->M());
151  }
152  }
153 
154  // ***
155  // *** do MASS CONSTRAINT FIT (J/psi)
156  // ***
157  for (j=0;j<jpsi.GetLength();++j)
158  {
159  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
160  mfitter.AddMassConstraint(m0_jpsi); // add the mass constraint
161  mfitter.Fit(); // do fit
162 
163  double chi2_m = mfitter.GetChi2(); // get chi2 of fit
164  double prob_m = mfitter.GetProb(); // access probability of fit
165  hjpsi_chi2_mf->Fill(chi2_m);
166  hjpsi_prob_mf->Fill(prob_m);
167 
168  if ( prob_m > 0.01 ) // when good enough, fill some histo
169  {
170  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
171  hjpsim_mcf->Fill(jfit->M());
172  }
173  }
174 
175  }
176 
177  // *** write out all the histos
178  out->cd();
179 
180  hjpsim_all->Write();
181  hpsim_all->Write();
182 
183  hjpsim_vf->Write();
184  hjpsim_4cf->Write();
185  hjpsim_mcf->Write();
186 
187  hjpsi_chi2_vf->Write();
188  hpsi_chi2_4c->Write();
189  hjpsi_chi2_mf->Write();
190 
191  hjpsi_prob_vf->Write();
192  hpsi_prob_4c->Write();
193  hjpsi_prob_mf->Write();
194 
195  hvpos->Write();
196 
197  out->Save();
198 
199 }
Int_t GetEntries()
Int_t i
Definition: run_full.C:25
TVector3 Pos() const
Definition: RhoCandidate.h:186
Int_t GetLength() const
Definition: RhoCandList.h:46
TH2F * hvpos
Definition: plot_eta_c.C:47
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
CandList piplus
FairRunAna * fRun
Definition: hit_dirc.C:58
CandList muplus
void Combine(RhoCandList &l1, RhoCandList &l2)
void Select(RhoParticleSelectorBase *pidmgr)
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
TFile * out
Definition: reco_muo.C:20
Double_t M() const
CandList piminus
Int_t GetEvent(Int_t n=-1)
CandList muminus