FairRoot/PandaRoot
ana_check_eta.C
Go to the documentation of this file.
1 class TCandList;
2 class TCandidate;
3 class TFitParams;
4 
5 int ana_check(int nevts=0)
6 {
7  TString OutFile="output.root";
8 
9  gStyle->SetOptFit(1011);
10  gROOT->Macro("$VMCWORKDIR/gconfig/rootlogon.C");
11 
12  TString inPidFile = "evt_pid.root";
13  TString inSimFile = "evt_points.root";
14 
15  TFile *inFile = TFile::Open(inSimFile,"READ");
16  TTree *tree=(TTree *) inFile->Get("pndsim") ;
17  tree->AddFriend("pndsim",inPidFile);
18 
19  TClonesArray* mc_array=new TClonesArray("PndMCTrack");
20  tree->SetBranchAddress("MCTrack",&mc_array);
21 
22  FairMCEventHeader* evthead;
23  tree->SetBranchAddress("MCEventHeader.", &evthead);
24 
25  TFile *out = TFile::Open(OutFile,"RECREATE");
26 
27  // the PndEventReader takes care about file/event handling
28  PndEventReader evr(inPidFile);
29 
30  TH1F *h_etac_nocut=new TH1F("h_etac_nocut","m(eta_c), (no cuts);E, GeV",100,2.5,3.5);
31  TH1F *h_mphi_nocuts=new TH1F("h_mphi_nocuts","#phi: m(K+ K-) (no cuts)",100,0.95,1.1);
32 
33  TH1F *h_etac_pid=new TH1F("h_etac_pid","m(eta_c), (MC PID);E, GeV",100,2.5,3.5);
34  TH1F *h_mphi_pid=new TH1F("h_mphi_pid","#phi: m(K+ K-) (MC PID)",100,0.95,1.1);
35 
36  TH1F *h_etac_4c=new TH1F("h_etac_4c","m(eta_c), 4C-fit",100,2.5,3.5);
37  TH1F *h_mphi_4c=new TH1F("h_mphi_4c","#phi: m(K+ K-) (4C-fit)",100,0.95,1.1);
38 
39  TH1F *h_etac_vtx=new TH1F("h_etac_vtx","m(eta_c), Vertex fit",100,2.5,3.5);
40  TH1F *h_mphi_vtx=new TH1F("h_mphi_vtx","#phi: m(K+ K-) (Vertex fit)",100,0.95,1.1);
41 
42  TH1F *h_etac_phimass=new TH1F("h_etac_phimass","m(eta_c), (cut on #phi mass);E, GeV",100,2.5,3.5);
43  TH1F *h_mphi_final=new TH1F("h_mphi_final","#phi: m(K+ K-)",100,0.95,1.1);
44 
45  TH1F *nc=new TH1F("nc","n charged",20,0,20);
46 
47  TH1F *h_chi2_4c=new TH1F("h_chi2_4c","#chi^{2} 4C-fit;#chi^{2}/N_{df}",100,0,100);
48  TH1F *h_chi2b_4c=new TH1F("h_chi2b_4c","#chi^{2} 4C-fit;#chi^{2}/N_{df}",100,0,100);
49  TH1F *h_chi2_vtx=new TH1F("h_chi2_vtx","#chi^{2} vertex;#chi^{2}/N_{df}",100,0,100);
50  TH2F *hvpos = new TH2F("hvpos","(x,y) projection of fitted decay vertex",100,-5,5,100,-5,5);
51  TH1F *hvzpos = new TH1F("hvzpos","z position of fitted decay vertex",100,-10,10);
52  TH1F *hvtxresX = new TH1F("hvtxresX","X resolution of fitted decay vertex",100,-0.1,0.1);
53  TH1F *hvtxresY = new TH1F("hvtxresY","Y resolution of fitted decay vertex",100,-0.1,0.1);
54  TH1F *hvtxresZ = new TH1F("hvtxresZ","Z resolution of fitted decay vertex",100,-0.1,0.1);
55 
56  TPidMassSelector *phiMassSel=new TPidMassSelector("phi",1.02,0.02);
57 
58  TPidPlusSelector *kplusSel=new TPidPlusSelector("kplus");
59  TPidMinusSelector *kminusSel=new TPidMinusSelector("kminus");
60 
61  // the candidates lists we need
62  TCandList p1, p2, phi1, phi1_pid, etac, etac_pid, etac_nocut;
63 
64  int n_reco=0;
65  // Number of events in file and number of reconstructed eta_c to store in root file
66  TH1F *n_events=new TH1F("n_events","total number of events",1,0,1);
67  TH1F *n_etac=new TH1F("n_etac","number of reconstructed eta_c",1,0,1);
68 
69  TLorentzVector ini(0,0,3.6772,4.7333);
70 
71  if (nevts==0) nevts=evr.GetEntries();
72 
73  n_events->SetBinContent(1,nevts);
74  // cout << "nevts " << nevts << "\n";
75  int i=0,j=0, k=0, l=0;
76 
77  // *************
78  // this is the loop through the events ... as simple as this...
79  // ****************
80  while (evr.GetEvent() && i++<nevts)
81  {
82  //cout << "evt: " << i << endl;
83  if ((i%100)==0) cout<<"evt " << i << "\n";
84  //if (!((i+1)%100)) cout<<"evt " << i << "\n";
85  evr.FillList(p1,"Charged");
86  evr.FillList(p2,"Charged");
87 
88  p1.Select(kplusSel);
89  p2.Select(kminusSel);
90 
91  int nchrg=p1.GetLength()+p2.GetLength();
92  nc->Fill(nchrg);
93 
94  for (j=0;j<p1.GetLength();++j) {
95  p1[j].SetMass(TRho::Instance()->GetPDG()->GetParticle(321)->Mass());
96  }
97  for (j=0;j<p2.GetLength();++j) {
98  p2[j].SetMass(TRho::Instance()->GetPDG()->GetParticle(321)->Mass());
99  }
100 
101  phi1.Combine(p1,p2);
102 
103  for (j=0;j<phi1.GetLength();++j) h_mphi_nocuts->Fill(phi1[j].M());
104 
105  phi1.Select(phiMassSel);
106  etac_nocut.Combine(phi1,phi1);
107 
108  for (l=0;l<etac_nocut.GetLength();++l) {
109  h_etac_nocut->Fill(etac_nocut[l].M());
110  }
111 
112  tree->GetEntry(i-1);
113  TVector3 mcVertex, mcD1Vertex, mcD2vertex;
114  evthead->GetVertex(mcVertex);
115 
116  // MC PID
117  // Leave only kaons in particle lists
118  int n_removed=0;
119  int ii=0;
120  for (l=0;l<p1.GetLength();++l) {
121  ii=l-n_removed;
122  if (p1[ii].GetMicroCandidate().GetMcIndex()>-1){
123  PndMCTrack *mcTrack = (PndMCTrack*)mc_array->At(p1[ii].GetMicroCandidate().GetMcIndex());
124  if (mcTrack!=0)
125  {
126  if ((mcTrack->GetPdgCode()!=321) || (mcTrack->GetMotherID()!=-1))
127  {
128  p1.Remove(p1[ii]);
129  n_removed++;
130  }
131  if (mcTrack==0) mcD1Vertex = mcTrack->GetStartVertex();
132  if (mcTrack==2) mcD2Vertex = mcTrack->GetStartVertex();
133  }
134  else
135  {
136  std::cout<<"stt h: " << p1[ii].GetMicroCandidate().GetSttHits() << std::endl;
137  std::cout<<"Kaon list 1, element "<<l<<" has no assosiated mcTRack"<<std::endl;
138  }
139  }
140  }
141 
142  n_removed=0;
143  ii=0;
144  for (l=0;l<p2.GetLength();++l) {
145  ii=l-n_removed;
146  if (p2[ii].GetMicroCandidate().GetMcIndex()>-1){
147  PndMCTrack *mcTrack = (PndMCTrack*)mc_array->At(p2[ii].GetMicroCandidate().GetMcIndex());
148  if (mcTrack!=0)
149  {
150  if ((mcTrack->GetPdgCode()!=-321) || (mcTrack->GetMotherID()!=-1))
151  {
152  p2.Remove(p2[ii]);
153  n_removed++;
154  }
155  if (mcTrack==1) mcD1Vertex = mcTrack->GetStartVertex();
156  if (mcTrack==3) mcD2Vertex = mcTrack->GetStartVertex();
157  }
158  else
159  {
160  std::cout<<"stt h: " << p2[ii].GetMicroCandidate().GetSttHits() << std::endl;
161  std::cout<<"Kaon list 2, element "<<l<<" has no assosiated mcTRack"<<std::endl;
162  }
163  }
164  }
165 
166  phi1_pid.Combine(p1,p2);
167 
168  for (j=0;j<phi1_pid.GetLength();++j) h_mphi_pid->Fill(phi1_pid[j].M());
169  phi1_pid.Select(phiMassSel);
170  etac.Combine(phi1_pid,phi1_pid);
171 
172  for (l=0;l<etac.GetLength();++l) {
173  h_etac_pid->Fill(etac[l].M());
174  }
175  cout << " ";// << endl;
176  //if (use4cfit)
177  {
179  int best_i=0;
180  double best_chi2=1000;
181  TCandidate *ccfit;// = new TCandidate();
182  double m_phi1, m_phi2;
183  for (l=0;l<etac.GetLength();++l) {
184  Rho4CFitter fitter(etac[l],ini);
185  fitter.FitConserveMasses();
186  double chi2=fitter.GetChi2();
187  if (chi2<best_chi2)
188  {
189  best_chi2=chi2;
190  best_i = l;
191  ccfit = (TCandidate*)fitter.FittedCand(etac[l]);
192  TCandidate *phi1best = fitter.FittedCand(*(etac[l].Daughter(0)));
193  TCandidate *phi2best = fitter.FittedCand(*(etac[l].Daughter(1)));
194  TCandidate *k1best = fitter.FittedCand(*(etac[l].Daughter(0)->Daughter(0)));
195  TCandidate *k2best = fitter.FittedCand(*(etac[l].Daughter(0)->Daughter(1)));
196  TCandidate *k3best = fitter.FittedCand(*(etac[l].Daughter(1)->Daughter(0)));
197  TCandidate *k4best = fitter.FittedCand(*(etac[l].Daughter(1)->Daughter(1)));
198  TLorentzVector tlvk1 = k1best->P4();
199  TLorentzVector tlvk2 = k2best->P4();
200  TLorentzVector tlvk3 = k3best->P4();
201  TLorentzVector tlvk4 = k4best->P4();
202  m_phi1= (tlvk1+tlvk2).M();
203  m_phi2= (tlvk3+tlvk4).M();
204  }
205  h_chi2_4c->Fill(chi2/9); // Ndf=3N-3=9
206  }
207 
208  if (etac.GetLength()!=0) cout << "evt: " << i << "\tbest: " << best_chi2 << "\tlen " << etac.GetLength() << endl;
209 
210  if(/*(best_chi2<270)&&*/(etac.GetLength()!=0))
211  {
212  h_chi2b_4c->Fill(best_chi2/9); // Ndf=3N-3=9
213  h_etac_4c->Fill(ccfit->M());
214  h_mphi_4c->Fill(m_phi1);
215  h_mphi_4c->Fill(m_phi2);
216  h_mphi_final->Fill(m_phi1);
217  h_mphi_final->Fill(m_phi2);
218  if (((m_phi1>1.02-0.02)&&(m_phi1<1.02+0.02))&&((m_phi2>1.02-0.02)&&(m_phi2<1.02+0.02)))
219  {
220  h_etac_phimass->Fill(etac[best_i].M());
221  if ((etac[best_i].M()>2.9)&&(etac[best_i].M()<3.06))
222  n_reco++;
223  }
224  }
225  }
226  // else // use vertex fit
227  {
228  TCandList etac_vtx;
229  TCandidate *k1, *k2, *k3, *k4, *phi1tmp, *phi2tmp, *etac_tmp;
230 
231  //Combine 4 kaons directly to candidates
232  for (j=0;j<etac.GetLength();++j)
233  {
234  phi1tmp=etac[j].Daughter(0);
235  phi2tmp=etac[j].Daughter(1);
236 
237  k1=phi1tmp->Daughter(0);
238  k2=phi1tmp->Daughter(1);
239  k3=phi2tmp->Daughter(0);
240  k4=phi2tmp->Daughter(1);
241 
242  etac_tmp=k1->Combine(*k2,*k3,*k4);
243  etac_vtx.Add(*etac_tmp);
244  }
245 
246  for (l=0;l<etac.GetLength();++l) {
247  h_etac_pid->Fill(etac_vtx[l].M());
248  }
249 
251  int best_i=0;
252  double best_chi2=1000;
253  TCandidate *etacfit_best=0;
254  TCandidate *k1fit_best, *k2fit_best, *k3fit_best, *k4fit_best;
255  TCandidate *phi1fit_best, *phi2fit_best;
256  TVector3 bestPos;
257  Float_t etacvtx_mass;
258  for (j=0;j<etac_vtx.GetLength();++j)
259  {
260  RhoKinVtxFitter vtxfitter(etac_vtx[j]); // instantiate a vertex fitter
261  vtxfitter.Fit(); // do the vertex fit
262 
263  TCandidate *etacfit=vtxfitter.FittedCand(etac_vtx[j]); // request the fitted EtaC candidate
264  TVector3 etacVtx=etacfit->Pos(); // and the decay vertex position
265  double chi2_vtx=vtxfitter.GlobalChi2();
266  h_chi2_vtx->Fill(chi2_vtx/5); // Number degree of freedom 2N-3=5
267  // plot mass and vtx x,y projection after fit
268  hvpos->Fill(etacVtx.X(),etacVtx.Y());
269  hvzpos->Fill(etacVtx.Z());
270  if(chi2_vtx<best_chi2)
271  {
272  best_chi2=chi2;
273  best_i=l;
274  etacfit_best=etacfit;
275  k1fit_best=vtxfitter.FittedCand(*(etacfit_best->Daughter(0)));
276  k2fit_best=vtxfitter.FittedCand(*(etacfit_best->Daughter(1)));
277  k3fit_best=vtxfitter.FittedCand(*(etacfit_best->Daughter(2)));
278  k4fit_best=vtxfitter.FittedCand(*(etacfit_best->Daughter(3)));
279  etacvtx_mass = etacfit_best->M();
280  bestPos = etacfit->Pos();
281  }
282 
283  }
284  if(/*(best_chi2<150)&&*/(etac.GetLength()!=0))
285  {
286  //h_etac_vtx->Fill(etacfit_best->M());
287  h_etac_vtx->Fill(etacvtx_mass);
288  phi1fit_best=k1fit_best->Combine(*k2fit_best);
289  phi2fit_best=k3fit_best->Combine(*k4fit_best);
290  double m_phi1=phi1fit_best->M();
291  double m_phi2=phi2fit_best->M();
292  h_mphi_vtx->Fill(m_phi1);
293  h_mphi_vtx->Fill(m_phi2);
294  h_mphi_final->Fill(m_phi1);
295  h_mphi_final->Fill(m_phi2);
296  hvtxresX->Fill(mcVertex.X()-bestPos.X());
297  hvtxresY->Fill(mcVertex.Y()-bestPos.Y());
298  hvtxresZ->Fill(mcVertex.Z()-bestPos.Z());
299 
300  if (((m_phi1>1.02-0.03)&&(m_phi1<1.02+0.03))&&((m_phi2>1.02-0.03)&&(m_phi2<1.02+0.03)))
301  {
302  h_etac_phimass->Fill(etacfit_best->M());
303  if ((etacfit_best->M()>2.9)&&(etacfit_best->M()<3.06))
304  n_reco++;
305  }
306  }
307 
308  }
309  }
310  std::cout<<"Number of reconstructed eta_c = "<<n_reco<<std::endl;
311  n_etac->SetBinContent(1,n_reco);
312 
313 
314  out->cd();
315  n_etac->Write();
316  n_events->Write();
317  h_etac_nocut->Write();
318  h_etac_pid->Write();
319  h_etac_phimass->Write();
320  h_etac_vtx->Write();
321  h_etac_4c->Write();
322 
323  h_mphi_nocuts->Write();
324  h_mphi_pid->Write();
325  h_mphi_vtx->Write();
326  h_mphi_4c->Write();
327  h_mphi_final->Write();
328 
329  nc->Write();
330 
331  h_chi2_4c->Write();
332  h_chi2_vtx->Write();
333  hvzpos->Write();
334  hvpos->Write();
335  hvtxresX->Write();
336  hvtxresY->Write();
337  hvtxresZ->Write();
338  out->Save();
339 
340  return 0;
341 }
342 
343 void removeCombinatoric(TCandList &etac_list)
344 {
345  TCandidate *phi1, *phi2, *k1, *k2, *k3, *k4;
346  bool isOverlap;
347 
348  int len1=etac_list.GetLength();
349  int ii=0;
350 
351  for (int i=0; i<len1; ++i)
352  {
353  phi1=etac_list[ii].Daughter(0);
354  phi2=etac_list[ii].Daughter(1);
355  k1=phi1->Daughter(0);
356  k2=phi1->Daughter(1);
357  k3=phi2->Daughter(0);
358  k4=phi2->Daughter(1);
359 
360  if ((k1->IsCloneOf(*k3))||(k2->IsCloneOf(*k4)))
361  {
362  etac_list.Remove(etac_list[ii]);
363  std::cout<<"Overlap found, i="<<i<<std::endl;
364  }
365  else
366  ii++;
367  }
368 
369 }
370 
371 
TH1F * h_etac_4c
Definition: plot_eta_c.C:24
TH1F * h_etac_pid
Definition: plot_eta_c.C:18
TH1F * hvtxresX
TH1F * h_mphi_nocuts
Definition: plot_eta_c.C:27
Int_t i
Definition: run_full.C:25
TVectorD * n_etac
Definition: plot_eta_c.C:13
TTree * tree
Definition: plot_dirc.C:12
TH1F * nc
Definition: plot_eta_c.C:38
Bool_t FitConserveMasses()
Definition: Rho4CFitter.cxx:64
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
TH2F * hvpos
Definition: plot_eta_c.C:47
TClonesArray * mc_array
Definition: anaLmdCluster.C:25
int ana_check(int nevts=0)
Definition: ana_check_eta.C:5
void removeCombinatoric(TCandList &etac_list)
TH1F * h_mphi_4c
Definition: plot_eta_c.C:33
TH1F * h_mphi_pid
Definition: plot_eta_c.C:29
TH1F * h_etac_phimass
Definition: plot_eta_c.C:20
TString inFile
Definition: hit_dirc.C:8
TH1F * h_chi2_vtx
Definition: plot_eta_c.C:43
TH1F * h_chi2_4c
Definition: plot_eta_c.C:41
TString inSimFile
TH1F * h_chi2b_4c
TH1F * hvzpos
Definition: plot_eta_c.C:45
TH1F * h_etac_nocut
Definition: plot_eta_c.C:16
TFile * out
Definition: reco_muo.C:20
TH1F * h_mphi_final
Definition: plot_eta_c.C:35
TVectorD * n_events
Definition: plot_eta_c.C:14
TPad * p2
Definition: hist-t7.C:117
TH1F * hvtxresZ
TPad * p1
Definition: hist-t7.C:116
double GetChi2() const
Definition: RhoFitterBase.h:48
TH1F * h_mphi_vtx
Definition: plot_eta_c.C:31
TVector3 GetStartVertex() const
Definition: PndMCTrack.h:76
Int_t GetMotherID() const
Definition: PndMCTrack.h:74
TH1F * h_etac_vtx
Definition: plot_eta_c.C:22
TH1F * hvtxresY