FairRoot/PandaRoot
rho/PndTutAnaTask.cxx
Go to the documentation of this file.
1 // ************************************************************************
2 //
3 // psi(2S) -> J/psi (-> mu+ mu-) pi+ pi- Analysis Example Task
4 //
5 // for the Rho Tutorial, see
6 // http://panda-wiki.gsi.de/cgi-bin/viewauth/Computing/PandaRootRhoTutorial
7 //
8 // K.Goetzen 7/2013
9 // ************************************************************************
10 
11 
12 // The header file
13 #include "PndTutAnaTask.h"
14 
15 // C++ headers
16 #include <string>
17 #include <iostream>
18 
19 // FAIR headers
20 #include "FairRootManager.h"
21 #include "FairRunAna.h"
22 #include "FairRuntimeDb.h"
23 #include "FairRun.h"
24 #include "FairRuntimeDb.h"
25 
26 // ROOT headers
27 #include "TClonesArray.h"
28 #include "TVector3.h"
29 #include "TH1F.h"
30 #include "TH2F.h"
31 
32 // RHO headers
33 #include "RhoCandidate.h"
34 #include "RhoHistogram/RhoTuple.h"
35 #include "RhoFactory.h"
37 
38 // Analysis headers
39 #include "PndAnalysis.h"
40 #include "Rho4CFitter.h"
41 #include "RhoKinVtxFitter.h"
42 #include "RhoKinFitter.h"
43 #include "RhoVtxPoca.h"
44 
45 
46 using std::cout;
47 using std::endl;
48 
49 
50 // ----- Default constructor -------------------------------------------
52  FairTask("Panda Tutorial Analysis Task") {
53 }
54 // -------------------------------------------------------------------------
55 
56 
57 // ----- Destructor ----------------------------------------------------
59 // -------------------------------------------------------------------------
60 
61 
62 // ----- Method to select true PID candidates
64 {
65  int removed = 0;
66 
67  for (int ii=l.GetLength()-1;ii>=0;--ii)
68  {
69  if ( !(ana->McTruthMatch(l[ii])) )
70  {
71  l.Remove(l[ii]);
72  removed++;
73  }
74  }
75 
76  return removed;
77 }
78 // -------------------------------------------------------------------------
79 
80 
81 // ----- Public method Init --------------------------------------------
82 InitStatus PndTutAnaTask::Init()
83 {
84  // initialize analysis object
85  fAnalysis = new PndAnalysis();
86 
87  // reset the event counter
88  fEvtCount = 0;
89 
90  // Mass selector for the jpsi cands
91  fJpsiMassSel=new RhoMassParticleSelector("jpsi",3.096,1.0);
92 
93  // create the histograms
94  hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
95  hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
96 
97  hjpsim_lpid = new TH1F("hjpsim_lpid","J/#psi mass (loose pid)",200,0,4.5);
98  hpsim_lpid = new TH1F("hpsim_lpid","#psi(2S) mass (loose pid)",200,0,5);
99 
100  hjpsim_tpid = new TH1F("hjpsim_tpid","J/#psi mass (tight pid)",200,0,4.5);
101  hpsim_tpid = new TH1F("hpsim_tpid","#psi(2S) mass (tight pid)",200,0,5);
102 
103  hjpsim_trpid = new TH1F("hjpsim_trpid","J/#psi mass (true pid)",200,0,4.5);
104  hpsim_trpid = new TH1F("hpsim_trpid","#psi(2S) mass (true pid)",200,0,5);
105 
106  hjpsim_ftm = new TH1F("hjpsim_ftm","J/#psi mass (full truth match)",200,0,4.5);
107  hpsim_ftm = new TH1F("hpsim_ftm","#psi(2S) mass (full truth match)",200,0,5);
108 
109  hjpsim_nm = new TH1F("hjpsim_nm","J/#psi mass (no truth match)",200,0,4.5);
110  hpsim_nm = new TH1F("hpsim_nm","#psi(2S) mass (no truth match)",200,0,5);
111 
112  hjpsim_diff = new TH1F("hjpsim_diff","J/#psi mass diff to truth",100,-2,2);
113  hpsim_diff = new TH1F("hpsim_diff","#psi(2S) mass diff to truth",100,-2,2);
114 
115  hjpsim_vf = new TH1F("hjpsim_vf","J/#psi mass (vertex fit)",200,0,4.5);
116  hjpsim_4cf = new TH1F("hjpsim_4cf","J/#psi mass (4C fit)",200,0,4.5);
117  hjpsim_mcf = new TH1F("hjpsim_mcf","J/#psi mass (mass constraint fit)",200,0,4.5);
118 
119  hjpsi_chi2_vf = new TH1F("hjpsi_chi2_vf", "J/#psi: #chi^{2} vertex fit",100,0,10);
120  hpsi_chi2_4c = new TH1F("hpsi_chi2_4c", "#psi(2S): #chi^{2} 4C fit",100,0,250);
121  hjpsi_chi2_mf = new TH1F("hjpsi_chi2_mf", "J/#psi: #chi^{2} mass fit",100,0,10);
122 
123  hjpsi_prob_vf = new TH1F("hjpsi_prob_vf", "J/#psi: Prob vertex fit",100,0,1);
124  hpsi_prob_4c = new TH1F("hpsi_prob_4c", "#psi(2S): Prob 4C fit",100,0,1);
125  hjpsi_prob_mf = new TH1F("hjpsi_prob_mf", "J/#psi: Prob mass fit",100,0,1);
126 
127  hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-2,2,100,-2,2);
128 
129 
130  // *** the lorentz vector of the initial psi(2S)
131  fIni.SetXYZT(0, 0, 6.231552, 7.240065);
132 
133  return kSUCCESS;
134 }
135 
136 // -------------------------------------------------------------------------
137 
139 {
140  // Get run and runtime database
141  FairRun* run = FairRun::Instance();
142  if ( ! run ) Fatal("SetParContainers", "No analysis run");
143 }
144 
145 // -------------------------------------------------------------------------
146 
147 
148 // ----- Public method Exec --------------------------------------------
149 void PndTutAnaTask::Exec(Option_t*)
150 {
151  // *** some variables
152  int i=0,j=0, k=0, l=0;
153 
154  // necessary to read the next event
156 
157  if (!(++fEvtCount%100)) cout << "evt "<<fEvtCount<<endl;
158 
159  // *** RhoCandLists for the analysis
160  RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
161 
162 
163  // *** Select with no PID info ('All'); type and mass are set
164  fAnalysis->FillList(muplus, "MuonAllPlus");
165  fAnalysis->FillList(muminus, "MuonAllMinus");
166  fAnalysis->FillList(piplus, "PionAllPlus");
167  fAnalysis->FillList(piminus, "PionAllMinus");
168 
169  // *** combinatorics for J/psi -> mu+ mu-
170  jpsi.Combine(muplus, muminus);
171 
172 
173  // ***
174  // *** do the TRUTH MATCH for jpsi
175  // ***
176  jpsi.SetType(443);
177 
178  for (j=0;j<jpsi.GetLength();++j)
179  {
180  hjpsim_all->Fill( jpsi[j]->M() );
181 
182  if (fAnalysis->McTruthMatch(jpsi[j]))
183  {
184  hjpsim_ftm->Fill( jpsi[j]->M() );
185  hjpsim_diff->Fill( jpsi[j]->GetMcTruth()->M() - jpsi[j]->M() );
186  }
187  else
188  hjpsim_nm->Fill( jpsi[j]->M() );
189  }
190 
191  // ***
192  // *** do VERTEX FIT (J/psi)
193  // ***
194  for (j=0;j<jpsi.GetLength();++j)
195  {
196  RhoKinVtxFitter vtxfitter(jpsi[j]); // instantiate a vertex fitter
197  vtxfitter.Fit();
198 
199  double chi2_vtx = vtxfitter.GetChi2(); // access chi2 of fit
200  double prob_vtx = vtxfitter.GetProb(); // access probability of fit
201  hjpsi_chi2_vf->Fill(chi2_vtx);
202  hjpsi_prob_vf->Fill(prob_vtx);
203 
204  if ( prob_vtx > 0.01 ) // when good enough, fill some histos
205  {
206  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
207  TVector3 jVtx=jfit->Pos(); // and the decay vertex position
208 
209  hjpsim_vf->Fill(jfit->M());
210  hvpos->Fill(jVtx.X(),jVtx.Y());
211  }
212  }
213 
214  // *** some rough mass selection
215  jpsi.Select(fJpsiMassSel);
216 
217  // *** combinatorics for psi(2S) -> J/psi pi+ pi-
218  psi2s.Combine(jpsi, piplus, piminus);
219 
220 
221  // ***
222  // *** do the TRUTH MATCH for psi(2S)
223  // ***
224  psi2s.SetType(88888);
225  for (j=0;j<psi2s.GetLength();++j)
226  {
227  hpsim_all->Fill( psi2s[j]->M() );
228 
229  if (fAnalysis->McTruthMatch(psi2s[j]))
230  {
231  hpsim_ftm->Fill( psi2s[j]->M() );
232  hpsim_diff->Fill( psi2s[j]->GetMcTruth()->M() - psi2s[j]->M() );
233  }
234  else
235  hpsim_nm->Fill( psi2s[j]->M() );
236  }
237 
238  // ***
239  // *** do 4C FIT (initial psi(2S) system)
240  // ***
241  for (j=0;j<psi2s.GetLength();++j)
242  {
243  RhoKinFitter fitter(psi2s[j]); // instantiate the kin fitter in psi(2S)
244  fitter.Add4MomConstraint(fIni); // set 4 constraint
245  fitter.Fit(); // do fit
246 
247  double chi2_4c = fitter.GetChi2(); // get chi2 of fit
248  double prob_4c = fitter.GetProb(); // access probability of fit
249  hpsi_chi2_4c->Fill(chi2_4c);
250  hpsi_prob_4c->Fill(prob_4c);
251 
252  if ( prob_4c > 0.01 ) // when good enough, fill some histo
253  {
254  RhoCandidate *jfit = psi2s[j]->Daughter(0)->GetFit(); // get fitted J/psi
255 
256  hjpsim_4cf->Fill(jfit->M());
257  }
258  }
259 
260 
261  // ***
262  // *** do MASS CONSTRAINT FIT (J/psi)
263  // ***
264  for (j=0;j<jpsi.GetLength();++j)
265  {
266  RhoKinFitter mfitter(jpsi[j]); // instantiate the RhoKinFitter in psi(2S)
267  mfitter.AddMassConstraint(3.0965); // add the mass constraint
268  mfitter.Fit(); // do fit
269 
270  double chi2_m = mfitter.GetChi2(); // get chi2 of fit
271  double prob_m = mfitter.GetProb(); // access probability of fit
272  hjpsi_chi2_mf->Fill(chi2_m);
273  hjpsi_prob_mf->Fill(prob_m);
274 
275  if ( prob_m > 0.01 ) // when good enough, fill some histo
276  {
277  RhoCandidate *jfit = jpsi[j]->GetFit(); // access the fitted cand
278  hjpsim_mcf->Fill(jfit->M());
279  }
280  }
281 
282 
283  // ***
284  // *** TRUE PID combinatorics
285  // ***
286 
287  // *** do MC truth match for PID type
288  SelectTruePid(fAnalysis, muplus);
289  SelectTruePid(fAnalysis, muminus);
290  SelectTruePid(fAnalysis, piplus);
291  SelectTruePid(fAnalysis, piminus);
292 
293  // *** all combinatorics again with true PID
294  jpsi.Combine(muplus, muminus);
295  for (j=0;j<jpsi.GetLength();++j) hjpsim_trpid->Fill( jpsi[j]->M() );
296  jpsi.Select(fJpsiMassSel);
297 
298  psi2s.Combine(jpsi, piplus, piminus);
299  for (j=0;j<psi2s.GetLength();++j) hpsim_trpid->Fill( psi2s[j]->M() );
300 
301 
302  // ***
303  // *** LOOSE PID combinatorics
304  // ***
305 
306  // *** and again with PidAlgoMvd;PidAlgoStt;PidAlgoDrc and loose selection
307  fAnalysis->FillList(muplus, "MuonLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
308  fAnalysis->FillList(muminus, "MuonLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
309  fAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
310  fAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
311 
312  jpsi.Combine(muplus, muminus);
313  for (j=0;j<jpsi.GetLength();++j) hjpsim_lpid->Fill( jpsi[j]->M() );
314  jpsi.Select(fJpsiMassSel);
315 
316  psi2s.Combine(jpsi, piplus, piminus);
317  for (j=0;j<psi2s.GetLength();++j) hpsim_lpid->Fill( psi2s[j]->M() );
318 
319 
320  // ***
321  // *** TIGHT PID combinatorics
322  // ***
323 
324  // *** and again with PidAlgoMvd;PidAlgoStt and tight selection
325  fAnalysis->FillList(muplus, "MuonTightPlus", "PidAlgoMdtHardCuts");
326  fAnalysis->FillList(muminus, "MuonTightMinus", "PidAlgoMdtHardCuts");
327  fAnalysis->FillList(piplus, "PionLoosePlus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
328  fAnalysis->FillList(piminus, "PionLooseMinus", "PidAlgoMvd;PidAlgoStt;PidAlgoDrc");
329 
330  jpsi.Combine(muplus, muminus);
331  for (j=0;j<jpsi.GetLength();++j) hjpsim_tpid->Fill( jpsi[j]->M() );
332  jpsi.Select(fJpsiMassSel);
333 
334  psi2s.Combine(jpsi, piplus, piminus);
335  for (j=0;j<psi2s.GetLength();++j) hpsim_tpid->Fill( psi2s[j]->M() );
336 
337 }
338 
339 
341 {
342  hjpsim_all->Write();
343  hpsim_all->Write();
344  hjpsim_lpid->Write();
345  hpsim_lpid->Write();
346  hjpsim_tpid->Write();
347  hpsim_tpid->Write();
348  hjpsim_trpid->Write();
349  hpsim_trpid->Write();
350 
351  hjpsim_ftm->Write();
352  hpsim_ftm->Write();
353  hjpsim_nm->Write();
354  hpsim_nm->Write();
355 
356  hpsim_diff->Write();
357  hjpsim_diff->Write();
358 
359  hjpsim_vf->Write();
360  hjpsim_4cf->Write();
361  hjpsim_mcf->Write();
362 
363  hjpsi_chi2_vf->Write();
364  hpsi_chi2_4c->Write();
365  hjpsi_chi2_mf->Write();
366 
367  hjpsi_prob_vf->Write();
368  hpsi_prob_4c->Write();
369  hjpsi_prob_mf->Write();
370 
371  hvpos->Write();
372 
373 }
374 
virtual void SetParContainers()
Int_t run
Definition: autocutx.C:47
Int_t i
Definition: run_full.C:25
TVector3 Pos() const
Definition: RhoCandidate.h:186
Int_t GetLength() const
Definition: RhoCandList.h:46
Bool_t FillList(RhoCandList &l, TString listkey="All", TString pidTcaNames="", int trackHypothesis=-1)
TLorentzVector fIni
CandList piplus
virtual void Exec(Option_t *opt)
CandList muplus
void GetEventInTask()
void Combine(RhoCandList &l1, RhoCandList &l2)
PndAnalysis * fAnalysis
void Select(RhoParticleSelectorBase *pidmgr)
void SetType(const TParticlePDG *pdt, int start=0)
Double_t M() const
Int_t Remove(RhoCandidate *)
virtual void Finish()
RhoMassParticleSelector * fJpsiMassSel
CandList piminus
ClassImp(PndAnaContFact)
Bool_t McTruthMatch(RhoCandidate *cand, Int_t level=2, bool verbose=false)
int SelectTruePid(PndAnalysis *ana, RhoCandList &l)
virtual InitStatus Init()
CandList muminus