]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG1/VZERO/extractThr.C
add missing include
[u/mrichter/AliRoot.git] / PWG1 / VZERO / extractThr.C
CommitLineData
c014f07d 1void extractThr(const char *filename, const char *fparname, Float_t centrValue = 10.0, Float_t semiCentrValue = 50.0)
2{
3 gROOT->LoadMacro(Form("%s/AliAnaVZEROTrigger.cxx+",
4 gSystem->pwd()));
5
6 AliAnaVZEROTrigger *task = new AliAnaVZEROTrigger;
7 task->Setup(fparname);
8
9 TFile *f = TFile::Open(filename);
10
11 TList *list = (TList*)f->Get("coutput");
12
13 TH2F *h1 = (TH2F*)list->FindObject("fV0Charge2d");
14 TH2F *h2 = (TH2F*)list->FindObject("fV0Charge2dPercent");
15 h2->Divide(h1);
16 TProfile *h1_pfx = h1->ProfileX("h1_pfx");
17 TProfile *h1_pfy = h1->ProfileY("h1_pfy");
18 new TCanvas;
19 gStyle->SetPalette(1,0);
20 h2->Draw("colz");
21 new TCanvas;
22 TF1 *f1 = new TF1("f1","[0]*x",0,task->GetCentCuts(0));
23 f1->SetParameter(0,task->GetRatio());
24 h1_pfx->Fit(f1,"WR+");
25 h1_pfx->Fit(f1,"WR+");
26 h1_pfx->Fit(f1,"WR+");
27 h1_pfx->Draw();
28
29 new TCanvas;
30 TF1 *f2 = new TF1("f2","x/[0]",0,task->GetCentCuts(0)*f1->GetParameter(0));
31 f2->SetParameter(0,1./task->GetRatio());
32 h1_pfy->Fit(f2,"WR+");
33 h1_pfy->Fit(f2,"WR+");
34 h1_pfy->Fit(f2,"WR+");
35 h1_pfy->Draw();
36
37 TH1F *h3 = (TH1F*)list->FindObject("fV0Percent");
38 new TCanvas;
39 h3->Draw();
b13c9bc7 40 Double_t nTotalEvts = (Double_t)h3->Integral(h3->FindBin(0.0),h3->FindBin(91.0))/0.9;
c014f07d 41
42 Int_t nb = task->GetNThr();
43 new TCanvas;
44 h3->Sumw2();
45 TH1F **hh = new TH1F*[nb];
46 TH1F **hhall = new TH1F*[nb];
47 TF1 **ff = new TF1*[nb];
48 Double_t *nent = new Double_t[nb];
b13c9bc7 49 Double_t *nBinEvts = new Double_t[nb];
c014f07d 50 Double_t *nentall = new Double_t[nb];
51 for(Int_t i = 0; i < nb; ++i) {
52 hh[i] = (TH1F*)list->FindObject(Form("fV0PercentBins_%d",i));
53 hhall[i] = (TH1F*)list->FindObject(Form("fV0PercentBinsAll_%d",i));
54
55 nent[i] = hh[i]->GetEntries();
b13c9bc7 56 nBinEvts[i] = (Double_t)hh[i]->Integral(hh[i]->FindBin(0.0),hh[i]->FindBin(91.0));
c014f07d 57 nentall[i] = hhall[i]->GetEntries();
58
59 hh[i]->Sumw2();
60 hh[i]->Divide(hh[i],h3,1,1,"B");
61 // TF1 *ff = new TF1(Form("ff_%d",i),"x > [0] ? TMath::Exp(-(x-[0])*(x-[0])/(2.*[1]*[1])) : 1.0",0,100);
62 ff[i] = new TF1(Form("ff_%d",i),"1.-1./(1.+TMath::Exp(-(x-[0])/[1]))",0,100);
63 ff[i]->SetLineColor(i%9+1);
64 printf("%d\n",i);
e5f60c5a 65 if (nent[i] > 100) {
c014f07d 66 Double_t par0 = hh[i]->GetBinCenter(hh[i]->FindLastBinAbove(0.9));
67 printf("%f\n",par0);
68 ff[i]->SetParameter(0,par0);
e5f60c5a 69 ff[i]->SetParameter(1,0.3);
c014f07d 70 hh[i]->Fit(ff[i],"R+");
71 hh[i]->Fit(ff[i],"R+");
72 hh[i]->Fit(ff[i],"R+");
73 }
74 }
75
76 new TCanvas;
77 Double_t *eff99 = new Double_t[nb];
78 Double_t *eff05 = new Double_t[nb];
79 Double_t *purity = new Double_t[nb];
80 Double_t *purity2 = new Double_t[nb];
81 Bool_t first = kTRUE;
82 for(Int_t i = 0; i < nb; ++i) {
83 hh[i]->SetLineColor(i%9+1);
84 hh[i]->SetMarkerColor(i%9+1);
85 hh[i]->SetLineWidth(2.0);
86 hh[i]->SetStats(0);
87 if (first) {
88 first = kFALSE;
89 hh[i]->Draw("e");
90 hh[i]->GetXaxis()->SetTitle("Centrality percentile");
91 hh[i]->GetYaxis()->SetTitle("Efficiency");
92 }
93 else
94 hh[i]->Draw("e same");
95 eff99[i] = ff[i]->GetX(0.99);
96 eff05[i] = ff[i]->GetX(0.05)-eff99[i];
97 purity[i] = ff[i]->Integral(0,eff99[i])/ff[i]->Integral(0,100);
98 purity2[i] = purity[i]*nent[i]/nentall[i];
99 // printf("%d %d %d\n",i,nent[i],nentall[i]);
100 }
101
102 new TCanvas;
103 TGraph *gr99 = new TGraph(nb,eff99,eff05);
104 gr99->SetMarkerStyle(21);
105 gr99->GetHistogram()->GetXaxis()->SetTitle("Centrality percentile @ Eff=99%");
106 gr99->GetHistogram()->GetYaxis()->SetTitle("#Delta(Centrality percentile @ Eff=5% - @ Eff=99%)");
107 gr99->GetHistogram()->GetYaxis()->SetRangeUser(-7,10);
108 gr99->SetTitle("");
109 gr99->Draw("AP");
110
111 TLegend *leg0 = new TLegend(0.75,0.55,0.89,0.89);
112 leg0->AddEntry(gr99,"A and C sides","P");
113 leg0->Draw();
114
115 new TCanvas;
116 TGraph *pur = new TGraph(nb,eff99,purity);
117 pur->SetMarkerStyle(21);
118 pur->GetHistogram()->GetXaxis()->SetTitle("Centrality percentile @ Eff=99%");
119 pur->GetHistogram()->GetYaxis()->SetTitle("Purity");
120 pur->SetTitle("");
121 TF1 *fitPur = new TF1("fitPur","[0]-[1]*TMath::Exp(-x/[2])",6,50);
122 fitPur->SetParameter(0,0.9);
123 fitPur->SetParameter(1,0.4);
124 fitPur->SetParameter(2,9.0);
125 pur->Fit(fitPur,"R+");
126 pur->Fit(fitPur,"R+");
127 pur->Fit(fitPur,"R+");
128 pur->Draw("AP");
129
130 TGraph *pur2 = new TGraph(nb,eff99,purity2);
131 pur2->SetMarkerStyle(22);
132 pur2->GetHistogram()->GetXaxis()->SetTitle("Centrality percentile @ Eff=99%");
133 pur2->GetHistogram()->GetYaxis()->SetTitle("Purity (no phys & centr selection)");
134 pur2->SetTitle("");
135 TF1 *fitPur2 = new TF1("fitPur2","[0]-[1]*TMath::Exp(-x/[2])",6,50);
136 fitPur2->SetParameter(0,0.9);
137 fitPur2->SetParameter(1,0.4);
138 fitPur2->SetParameter(2,9.0);
139 pur2->Fit(fitPur2,"R+");
140 pur2->Fit(fitPur2,"R+");
141 pur2->Fit(fitPur2,"R+");
142 pur2->Draw("Psame");
143
144 TLegend *leg = new TLegend(0.75,0.55,0.89,0.89);
145 leg->AddEntry(pur,"A and C, Phys & centr selection","P");
146 leg->AddEntry(pur2,"A and C, No event selection","P");
147 leg->Draw();
148
149 new TCanvas;
150 Double_t *thrA = new Double_t[nb];
151 Double_t *thrC = new Double_t[nb];
152 Double_t ratio = task->GetRatio();
153 for(Int_t i = 0; i < nb; ++i) {
154 thrA[i] = task->GetThrA(i);
155 thrC[i] = task->GetThrC(i);
156 }
157 TGraph *thr = new TGraph(nb,eff99,thrA);
158 thr->SetMarkerStyle(21);
159 TF1 *fThr = new TF1("fThr","[0]*x*x*x+[1]*x*x+[2]*x+[3]",1.5,50);
160 thr->Fit(fThr,"R");
161 thr->Draw("AP");
162
b13c9bc7 163 new TCanvas;
164 Double_t *intThr = new Double_t[nb];
165 for(Int_t i = 0; i < nb; ++i) {
166 intThr[i] = 100.*nBinEvts[i]/nTotalEvts;
167 }
168 TGraph *grInt = new TGraph(nb,intThr,thrA);
169 grInt->SetMarkerStyle(21);
170 TF1 *fIntThr = new TF1("fIntThr","[0]*x*x*x+[1]*x*x+[2]*x+[3]",1.5,50);
171 grInt->Fit(fIntThr,"R");
172 grInt->Draw("AP");
173
e5f60c5a 174 new TCanvas;
c014f07d 175 TH2F *hall = (TH2F*)list->FindObject("fV0Charge2dAll");
176 h2->Draw("colz");
177 TLine *line0 = new TLine(fThr->Eval(centrValue),fThr->Eval(centrValue)*task->GetRatio(),25000,fThr->Eval(centrValue)*task->GetRatio());
178 line0->Draw("same");
179 TLine *line1 = new TLine(fThr->Eval(centrValue),fThr->Eval(centrValue)*task->GetRatio(),fThr->Eval(centrValue),25000);
180 line1->Draw("same");
181
b13c9bc7 182 TLine *line01 = new TLine(fIntThr->Eval(centrValue),fIntThr->Eval(centrValue)*task->GetRatio(),25000,fIntThr->Eval(centrValue)*task->GetRatio());
183 line01->Draw("same");
184 TLine *line11 = new TLine(fIntThr->Eval(centrValue),fIntThr->Eval(centrValue)*task->GetRatio(),fIntThr->Eval(centrValue),25000);
185 line11->Draw("same");
186
c014f07d 187 TLine *line2 = new TLine(fThr->Eval(semiCentrValue),fThr->Eval(semiCentrValue)*task->GetRatio(),25000,fThr->Eval(semiCentrValue)*task->GetRatio());
188 line2->Draw("same");
189 TLine *line3 = new TLine(fThr->Eval(semiCentrValue),fThr->Eval(semiCentrValue)*task->GetRatio(),fThr->Eval(semiCentrValue),25000);
190 line3->Draw("same");
191
192 printf("Slope param: %f %f\n",f1->GetParameter(0),f2->GetParameter(0));
193 printf("Mean slope param: %f\n",0.5*(f1->GetParameter(0)+f2->GetParameter(0)));
194 printf("Delta slope param: %f %%\n",100.*(f1->GetParameter(0)-f2->GetParameter(0))/(0.5*(f1->GetParameter(0)+f2->GetParameter(0))));
195
196 printf("A&C @ %.1f %% ThrA = %f ThrC = %f\n",centrValue,fThr->Eval(centrValue),fThr->Eval(centrValue)*task->GetRatio());
197 printf("A&C @ %.1f %% ThrA = %f ThrC = %f\n",semiCentrValue,fThr->Eval(semiCentrValue),fThr->Eval(semiCentrValue)*task->GetRatio());
198
b13c9bc7 199 printf("A&C (Integral) @ %.1f %% ThrA = %f ThrC = %f\n",centrValue,fIntThr->Eval(centrValue),fIntThr->Eval(centrValue)*task->GetRatio());
200 printf("A&C (Intergal) @ %.1f %% ThrA = %f ThrC = %f\n",semiCentrValue,fIntThr->Eval(semiCentrValue),fIntThr->Eval(semiCentrValue)*task->GetRatio());
201
c014f07d 202 printf("\n\n");
203 FILE *fout=fopen("./trigger.txt","w");
204 if (!fout) {
205 printf("Failed to open local result file\n");
206 return;
207 }
208 printf("%.0f %.0f %f %d %.0f %.0f %.0f %.0f\n\n",
209 task->GetMinThr(),
210 task->GetMaxThr(),
211 0.5*(f1->GetParameter(0)+f2->GetParameter(0)),
212 task->GetNThr(),
213 fThr->Eval(semiCentrValue),fThr->Eval(semiCentrValue)*task->GetRatio(),
b13c9bc7 214 fIntThr->Eval(centrValue),fIntThr->Eval(centrValue)*task->GetRatio());
c014f07d 215 fprintf(fout,"%.0f %.0f %f %d %.0f %.0f %.0f %.0f\n",
216 task->GetMinThr(),
217 task->GetMaxThr(),
218 0.5*(f1->GetParameter(0)+f2->GetParameter(0)),
219 task->GetNThr(),
220 fThr->Eval(semiCentrValue),fThr->Eval(semiCentrValue)*task->GetRatio(),
b13c9bc7 221 fIntThr->Eval(centrValue),fIntThr->Eval(centrValue)*task->GetRatio());
c014f07d 222 fclose(fout);
223}