]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFanalyzeMatching.C
Remove a run from LHC10e wo AODs; added run list for the new PbPb MC LHC11a10b_bis
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeMatching.C
CommitLineData
b213b8bd 1void AliTOFanalyzeMatching(const char* datafile)
2{
3
4 //
5 // Matching efficiency and contamination
6 // for different particle species
7 // (pions, kaons and protons).
8 // All histos are saved into a separate file
9 // datafile is assumed to be the file name containing
10 // the results of the matching in TNtuple format.
11 //
12 // Author: F. Pierella | pierella@bo.infn.it
13 //
14 // Use case:
15 // start root
16 // root[0] .L AliTOFanalyzeMatching.C
17 // root[1] AliTOFanalyzeMatching("matchingNtuple.root")
18
19 // output (histos!) filename
20 char outFileName[100];
21 strcpy(outFileName,"histo");
22 strcat(outFileName,datafile);
23
24 // dummy histos (for normalization)
25 TH1F* hpitot= new TH1F("hpitot","",12,0.,3.);
26 TH1F* hkatot= new TH1F("hkatot","",12,0.,3.);
27 TH1F* hprtot= new TH1F("hprtot","",12,0.,3.);
28
29 TH1F* hpimatched= new TH1F("hpimatched","",12,0.,3.);
30 TH1F* hkamatched= new TH1F("hkamatched","",12,0.,3.);
31 TH1F* hprmatched= new TH1F("hprmatched","",12,0.,3.);
32
33
34 // matching efficiency histos
35 TH1F* hpimatcheff= new TH1F("hpimatcheff","Matching efficiency for pions",12,0.,3.);
36 TH1F* hkamatcheff= new TH1F("hkamatcheff","Matching efficiency for kaons",12,0.,3.);
37 TH1F* hprmatcheff= new TH1F("hprmatcheff","Matching efficiency for protons",12,0.,3.);
38
39 // matching contamination histos
40 TH1F* hpimatchcon= new TH1F("hpimatchcon","Matching contamination for pions",12,0.,3.);
41 TH1F* hkamatchcon= new TH1F("hkamatchcon","Matching contamination for kaons",12,0.,3.);
42 TH1F* hprmatchcon= new TH1F("hprmatchcon","Matching contamination for protons",12,0.,3.);
43
44
45 TFile *file = TFile::Open(datafile,"old");
46 TNtuple* fNtuple= (TNtuple*)file->Get("Ntuple"); // get ntuple from file
47 Int_t nvar = fNtuple->GetNvar(); cout <<"N of var.="<< nvar << endl;
48 fNtuple->GetEvent(0);
49
50 file->cd();
51 Int_t nparticles = (Int_t)fNtuple->GetEntries();
52
53 for (Int_t i=0; i < nparticles; i++) {
54 fNtuple->GetEvent(i);
55 Int_t event=fNtuple->GetLeaf("event")->GetValue();
56 Int_t pdgcode=fNtuple->GetLeaf("ipart")->GetValue();
57 Int_t matc=fNtuple->GetLeaf("matc")->GetValue(0);
58 Float_t px=fNtuple->GetLeaf("pxvtx")->GetValue(0);
59 Float_t py=fNtuple->GetLeaf("pyvtx")->GetValue(0);
60 Float_t pz=fNtuple->GetLeaf("pzvtx")->GetValue(0);
61
62 Float_t pvtx=TMath::Sqrt(px*px+py*py+pz*pz);
63 Float_t ptvtx=TMath::Sqrt(px*px+py*py);
64 Int_t abspdgcode=TMath::Abs(pdgcode);
65
66
67 // N (1+2+3+4+(-4)) cases
68 if(matc>=1 || matc==-4){
69 switch(abspdgcode){
70 case 211:
71 hpitot->Fill(pvtx);
72 break;
73 case 321:
74 hkatot->Fill(pvtx);
75 break;
76 case 2212:
77 hprtot->Fill(pvtx);
78 break;
79 }
80 }
81
82
83 // N_matched (3+4) cases
84 if(matc==3 || matc==4){
85 switch(abspdgcode){
86 case 211:
87 hpimatched->Fill(pvtx);
88 break;
89 case 321:
90 hkamatched->Fill(pvtx);
91 break;
92 case 2212:
93 hprmatched->Fill(pvtx);
94 break;
95 }
96 }
97
98
99 // N_t (3) case
100 if(matc==3){
101 switch(abspdgcode){
102 case 211:
103 hpimatcheff->Fill(pvtx);
104 break;
105 case 321:
106 hkamatcheff->Fill(pvtx);
107 break;
108 case 2212:
109 hprmatcheff->Fill(pvtx);
110 break;
111 }
112 }
113
114 // N_w (4) case
115 if(matc==4){
116 switch(abspdgcode){
117 case 211:
118 hpimatchcon->Fill(pvtx);
119 break;
120 case 321:
121 hkamatchcon->Fill(pvtx);
122 break;
123 case 2212:
124 hprmatchcon->Fill(pvtx);
125 break;
126 }
127 }
128
129 }
130
131 // histo normalization
132 // efficiency
133 hpimatcheff->Divide(hpitot);
134 hkamatcheff->Divide(hkatot);
135 hprmatcheff->Divide(hprtot);
136
137 // contamination
138 hpimatchcon->Divide(hpimatched);
139 hkamatchcon->Divide(hkamatched);
140 hprmatchcon->Divide(hprmatched);
141
142
143 TFile *houtfile = new TFile(outFileName,"recreate");
144 houtfile->cd();
145
146 hpitot->Write();
147 hkatot->Write();
148 hprtot->Write();
149
150 hpimatched->Write();
151 hkamatched->Write();
152 hprmatched->Write();
153
154 hpimatcheff->Write();
155 hkamatcheff->Write();
156 hprmatcheff->Write();
157
158 hpimatchcon->Write();
159 hkamatchcon->Write();
160 hprmatchcon->Write();
161 houtfile->Close();
162
163 cout << "File " << outFileName << " with histos has been created" << endl;
164}