]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGGA/PHOSTasks/PHOS_PbPb/macros/MakeMmixPi0.C
force directory creation
[u/mrichter/AliRoot.git] / PWGGA / PHOSTasks / PHOS_PbPb / macros / MakeMmixPi0.C
CommitLineData
b8bea077 1/* $Id$ */
2
8e6cc59c 3const Double_t kMean=0.136 ; //Approximate peak position to facilitate error estimate
4
5//-----------------------------------------------------------------------------
91b112bd 6void MakeMmixPi0(const TString filename = "Pi0Flow_000167920.root",
7 const Int_t centrality=0,
8 const char* pid="CPV")
8e6cc59c 9{
10 //Fit Real/Mixed ratio, normalize Mixed and subtract it from Real
11
91b112bd 12 TFile * f = new TFile(filename) ;
8b4ce769 13 //TList *histoList = (TList*)f->Get("PHOSPi0Flow");
14 TList *histoList = (TList*)f->Get("PHOSPi0Flow/PHOSPi0FlowCoutput1"); // lego train
15
8e6cc59c 16 char key[125] ;
17
18 TH1F * hev = (TH1F*)histoList->FindObject("hTotSelEvents") ;
19 TH2F * hCentrality = (TH2F*)histoList->FindObject("hCenPHOSCells") ;
20 TH1D * hCentrality1 = hCentrality->ProjectionX();
21
22 printf("TotSelEvents: %f \n",hev->GetBinContent(4)) ;
23 printf("Centrality: %f \n",hCentrality1->Integral()) ;
24
25 TString inputKey;
26 TString outputKey = Form("%s_cent%d",pid,centrality);
27
28 TH2F *h , *hAdd;
29 TH2F *hm, *hmAdd;
30
31 if (centrality >= 0 && centrality < 4) {
32 printf("\tCentrality %d\n",centrality);
33 inputKey = Form("hPi0%s_cen%d" ,pid,centrality);
34 TH2F *h = (TH2F*)histoList->FindObject(inputKey) ;
35 inputKey = Form("hMiPi0%s_cen%d",pid,centrality);
36 TH2F *hm= (TH2F*)histoList->FindObject(inputKey) ;
37 if (h==0) {
38 printf("Missing histogram %s\n",inputKey);
39 return;
40 }
41 }
42 else {
43 printf("Wrong centrality %d. Allowed values are 0,1,2,3,10.\n",centrality);
44 return;
45 }
46
91b112bd 47 Int_t nPadX = 3, nPadY = 2;
48 Int_t nbin=6 ;
49 Double_t xa[21]={1., 2., 3., 4., 5., 7., 10.} ;
8e6cc59c 50
51 PPRstyle();
52 gStyle->SetPadLeftMargin(0.14);
53 gStyle->SetPadRightMargin(0.01);
54 gStyle->SetPadTopMargin(0.01);
55 gStyle->SetPadBottomMargin(0.08);
56
57 //Fit real only
58 //Linear Bg
59 char kkey[55];
60 sprintf(kkey,outputKey.Data()) ;
61 char key2[155];
62 sprintf(key,"Mix%s",kkey) ;
63 sprintf(key2,"%s_mr1",key) ;
64 TH1D * mr1 = new TH1D(key2,"Mass",nbin,xa) ;
65 sprintf(key2,"%s_sr1",key) ;
66 TH1D * sr1 = new TH1D(key2,"Width",nbin,xa) ;
67 sprintf(key2,"%s_ar1",key) ;
68 TH1D * ar1 = new TH1D(key2,"a",nbin,xa) ;
69 sprintf(key2,"%s_br1",key) ;
70 TH1D * br1 = new TH1D(key2,"a",nbin,xa) ;
71 sprintf(key2,"%s_yr1",key) ;
72 TH1D * nr1 = new TH1D(key2,"Raw yield",nbin,xa) ;
73 sprintf(key2,"%s_yr1int",key) ;
74 TH1D * nr1int = new TH1D(key2,"Raw yield, integrated",nbin,xa) ;
75
76 //Quadratic Bg
77 sprintf(key2,"%s_mr2",key) ;
78 TH1D * mr2 = new TH1D(key2,"Mass",nbin,xa) ;
79 sprintf(key2,"%s_sr2",key) ;
80 TH1D * sr2 = new TH1D(key2,"Width",nbin,xa) ;
81 sprintf(key2,"%s_ar2",key) ;
82 TH1D * ar2 = new TH1D(key2,"a",nbin,xa) ;
83 sprintf(key2,"%s_br2",key) ;
84 TH1D * br2 = new TH1D(key2,"a",nbin,xa) ;
85 sprintf(key2,"%s_cr2",key) ;
86 TH1D * cr2 = new TH1D(key2,"a",nbin,xa) ;
87 sprintf(key2,"%s_yr2",key) ;
88 TH1D * nr2 = new TH1D(key2,"Raw yield",nbin,xa) ;
89 sprintf(key2,"%s_yr2int",key) ;
90 TH1D * nr2int = new TH1D(key2,"Raw yield, integrated",nbin,xa) ;
91
92 TF1 * fit1 = new TF1("fit",CB,0.,1.,6) ;
93 fit1->SetParName(0,"A") ;
94 fit1->SetParName(1,"m_{0}") ;
95 fit1->SetParName(2,"#sigma") ;
96 fit1->SetParName(3,"a_{0}") ;
97 fit1->SetParName(4,"a_{1}") ;
98 fit1->SetParName(5,"a_{2}") ;
99 fit1->SetLineWidth(2) ;
100 fit1->SetLineColor(2) ;
101 TF1 * fgs = new TF1("gs",CBs,0.,1.,4) ;
102 fgs->SetParName(0,"A") ;
103 fgs->SetParName(1,"m_{0}") ;
104 fgs->SetParName(2,"#sigma") ;
105 fgs->SetParName(3,"B") ;
106 fgs->SetLineColor(2) ;
107 fgs->SetLineWidth(1) ;
108
109 TF1 * fit2 = new TF1("fit2",CB2,0.,1.,7) ;
110 fit2->SetParName(0,"A") ;
111 fit2->SetParName(1,"m_{0}") ;
112 fit2->SetParName(2,"#sigma") ;
113 fit2->SetParName(3,"a_{0}") ;
114 fit2->SetParName(4,"a_{1}") ;
115 fit2->SetParName(5,"a_{2}") ;
116 fit2->SetParName(6,"a_{3}") ;
117 fit2->SetLineWidth(2) ;
118 fit2->SetLineColor(4) ;
119 fit2->SetLineStyle(2) ;
120
121 TF1 * fbg1 = new TF1("bg1",BG1,0.,1.,3) ;
122 TF1 * fbg2 = new TF1("bg2",BG2,0.,1.,4) ;
123
124 TCanvas * c3 = new TCanvas("mggFit1_Signal","mggFitCB",10,10,1200,800) ;
125 c3->Divide(nPadX,nPadY) ;
126
127 TCanvas * c1 = new TCanvas("mggFit1","mggFit1",10,10,1200,800) ;
128 c1->Divide(nPadX,nPadY) ;
129 c1->cd(0) ;
130
131 TAxis * pta=h->GetYaxis() ;
132 TAxis * ma=h->GetXaxis() ;
133 for(Int_t i=1;i<=nbin;i++){
134 c1->cd(i) ;
135 printf("\n\t%.1f<pt<%.1f GeV/c\n",xa[i-1],xa[i]);
136 Int_t imin=pta->FindBin(xa[i-1]+0.0001);
137 Int_t imax=pta->FindBin(xa[i]-0.0001) ;
138 Double_t pt=(xa[i]+xa[i-1])/2. ;
139 TH1D * hp = h->ProjectionX("re",imin,imax) ;
140 hp->Sumw2() ;
141 TH1D * hpm= hm->ProjectionX("mi",imin,imax) ;
142 hpm->Sumw2() ;
143 if(i<=17){
144 hp ->Rebin(4) ;
145 hpm->Rebin(4) ;
146 }
147 else{
148 hp ->Rebin(5) ;
149 hpm->Rebin(5) ;
150 }
151 for(Int_t ib=1; ib<=hp->GetNbinsX();ib++){if(hp ->GetBinContent(ib)==0)hp ->SetBinError(ib,1.);}
152 for(Int_t ib=1; ib<=hp->GetNbinsX();ib++){if(hpm->GetBinContent(ib)==0)hpm->SetBinError(ib,1.);}
153 TH1D * hpm2 = (TH1D*)hpm->Clone("Bg1") ;
154 TH1D * hpcopy = (TH1D*)hp ->Clone("hpcopy") ;
155 TH1D * hp2 = (TH1D*)hp ->Clone("hp2") ;
156 hpcopy->SetXTitle("M_{#gamma#gamma} (GeV/c^{2})");
157 hp2 ->SetXTitle("M_{#gamma#gamma} (GeV/c^{2})");
158 hpcopy->Divide(hpm) ;
159 sprintf(key,"PID=%s, %3.1f<p_{T}<%3.1f GeV/c",pid,xa[i-1],xa[i]) ;
160 hpcopy->SetTitle(key) ;
161 hpcopy->SetMarkerStyle(20) ;
162 hpcopy->SetMarkerSize(0.7) ;
163
164 fit1->SetParameters(0.001,0.136,0.0055,0.0002,-0.002,0.0) ;
165 fit1->SetParLimits(0,0.000,1.000) ;
166 fit1->SetParLimits(1,0.120,0.145) ;
167 fit1->SetParLimits(2,0.005,0.012) ;
168
169 Double_t rangeMin=0.05 ;
170 Double_t rangeMax=0.3 ;
171 if(centrality==0) rangeMax=0.4 ;
172 if(i==1){
173 rangeMin=0.06 ;
174 rangeMax=0.25 ;
175 }
176 hpcopy->Fit(fit1,"Q" ,"",rangeMin,rangeMax) ;
177 hpcopy->Fit(fit1,"MQ","",rangeMin,rangeMax) ;
178
179 ar1->SetBinContent(i,fit1->GetParameter(3)) ;
180 ar1->SetBinError (i,fit1->GetParError(3)) ;
181 br1->SetBinContent(i,fit1->GetParameter(4)) ;
182 br1->SetBinError (i,fit1->GetParError(4)) ;
183
184 fit2->SetParameters(fit1->GetParameters()) ;
185 fit2->SetParLimits(0,0.000,1.000) ;
186 fit2->SetParLimits(1,0.120,0.145) ;
187 fit2->SetParLimits(2,0.005,0.012) ;
188
189 hpcopy->Fit(fit2,"+NQ","",rangeMin,rangeMax) ;
190 hpcopy->Fit(fit2,"+MQ","",rangeMin,rangeMax) ;
191
192 ar2->SetBinContent(i,fit2->GetParameter(3)) ;
193 ar2->SetBinError (i,fit2->GetParError(3)) ;
194 br2->SetBinContent(i,fit2->GetParameter(4)) ;
195 br2->SetBinError (i,fit2->GetParError(4)) ;
196 cr2->SetBinContent(i,fit2->GetParameter(5)) ;
197 cr2->SetBinError (i,fit2->GetParError(5)) ;
198 hpcopy->GetXaxis()->SetRangeUser(0.05,0.30) ;
199 hpcopy->Draw() ;
200 c1->Update() ;
201
202 fbg1->SetParameters(fit1->GetParameter(3),fit1->GetParameter(4),fit1->GetParameter(5));
203 fbg2->SetParameters(fit2->GetParameter(3),fit2->GetParameter(4),fit2->GetParameter(5),
204 fit2->GetParameter(6));
205
206 Double_t intRangeMin = PeakPosition(pt)-3.*PeakWidth(pt) ;
207 Double_t intRangeMax = PeakPosition(pt)+3.*PeakWidth(pt) ;
208 Int_t intBinMin = hp->GetXaxis()->FindBin(intRangeMin) ;
209 Int_t intBinMax = hp->GetXaxis()->FindBin(intRangeMax) ;
210 Double_t errStat = hpm->Integral(intBinMin,intBinMax);
211
212 hpm ->Multiply(fbg1) ;
213 hpm2->Multiply(fbg2) ;
214 hp ->Add(hpm ,-1.) ;
215 hp2 ->Add(hpm2,-1.) ;
216
217 c3->cd(i) ;
218
219 Int_t binPi0 = hp->FindBin(kMean);
220 fgs->SetParameters(hp->Integral(binPi0-1,binPi0+1)/3.,fit1->GetParameter(1),fit1->GetParameter(2)) ;
221 fgs->SetParLimits(0,0.000,1.e+5) ;
222 fgs->SetParLimits(1,0.120,0.145) ;
223 fgs->SetParLimits(2,0.004,0.02) ;
224 hp->Fit(fgs,"Q","",rangeMin,rangeMax) ;
225 hp->SetMaximum(hp2->GetMaximum()*1.5) ;
226 hp->SetMinimum(hp2->GetMinimum()*1.1) ;
227 hp->SetMarkerStyle(20) ;
228 hp->SetMarkerSize(0.7) ;
229 mr1->SetBinContent(i,fgs->GetParameter(1)) ;
230 mr1->SetBinError (i,fgs->GetParError(1) ) ;
231 sr1->SetBinContent(i,TMath::Abs(fgs->GetParameter(2))) ;
232 sr1->SetBinError (i,fgs->GetParError(2) ) ;
233
234 Double_t y=fgs->GetParameter(0)/hp->GetXaxis()->GetBinWidth(1) ;
235 nr1->SetBinContent(i,y) ;
236 Double_t ey=fgs->GetParError(0)/hp->GetXaxis()->GetBinWidth(1) ;
237 nr1->SetBinError(i,ey) ;
238
239 Double_t npiInt = hp->Integral(intBinMin,intBinMax) ;
240 Double_t norm = fbg1->GetParameter(0) ;
241 Double_t normErr= fbg1->GetParError(0) ;
242 if(npiInt>0.){
243 nr1int->SetBinContent(i,npiInt) ;
244 nr1int->SetBinError(i,TMath::Sqrt(npiInt + norm*errStat + normErr*normErr*errStat*errStat + norm*norm*errStat)) ;
245 }
246 hp2->GetXaxis()->SetRangeUser(0.05,0.3) ;
247 hp2->SetMaximum(hp2->GetMaximum()*1.5) ;
248 hp2->SetMinimum(hp2->GetMinimum()*1.1) ;
249 hp2->SetMarkerStyle(20) ;
250 hp2->SetMarkerSize(0.7) ;
251 hp2->Fit(fgs,"Q","",rangeMin,rangeMax) ;
252 mr2->SetBinContent(i,fgs->GetParameter(1)) ;
253 mr2->SetBinError (i,fgs->GetParError(1)) ;
254 sr2->SetBinContent(i,TMath::Abs(fgs->GetParameter(2))) ;
255 sr2->SetBinError (i,fgs->GetParError(2)) ;
256 y=fgs->GetParameter(0)/hp->GetXaxis()->GetBinWidth(1) ;
257 nr2->SetBinContent(i,y) ;
258 ey= fgs->GetParError(0)/hp->GetXaxis()->GetBinWidth(1) ;
259 nr2->SetBinError(i,ey) ;
260 npiInt=hp2->Integral(intBinMin,intBinMax) ;
261 norm=fbg2->GetParameter(0) ;
262 normErr=fbg2->GetParError(0) ;
263 if(npiInt>0.){
264 nr2int->SetBinContent(i,npiInt) ;
265 nr2int->SetBinError(i,TMath::Sqrt(npiInt + norm*errStat + normErr*normErr*errStat*errStat + norm*norm*errStat)) ;
266 }
267 hp2->SetTitle(key) ;
268 hp2->Draw() ;
269 c3->Update() ;
270
271 delete hp ;
272 delete hpm ;
273 delete hpm2 ;
274 }
275 char name[55] ;
276
277 sprintf(name,"Pi0_ratio_%s.eps",kkey) ;
278 c1->Print(name) ;
279 sprintf(name,"Pi0_signal_%s.eps",kkey) ;
280 c3->Print(name) ;
281
282 //Normalize by the number of events
283 Int_t cMin,cMax;
284 if (centrality == 0) {
285 cMin=1;
286 cMax=10;
287 }
288 else if (centrality == 1) {
289 cMin=21;
290 cMax=40;
291 }
292 else if (centrality == 2) {
293 cMin=41;
294 cMax=60;
295 }
296 else if (centrality == 3) {
297 cMin=61;
298 cMax=80;
299 }
300 Double_t nevents = hCentrality1->Integral(cMin,cMax);
301 nr1 ->Scale(1./nevents) ;
302 nr1int->Scale(1./nevents) ;
303 nr2 ->Scale(1./nevents) ;
304 nr2int->Scale(1./nevents) ;
305 printf("Nevents=%f \n",nevents) ;
306
307
308 TFile fout("LHC10h_Pi0_FitResult.root","update");
309 mr1->Write(0,TObject::kOverwrite) ;
310 sr1->Write(0,TObject::kOverwrite) ;
311 ar1->Write(0,TObject::kOverwrite) ;
312 br1->Write(0,TObject::kOverwrite) ;
313 nr1->Write(0,TObject::kOverwrite) ;
314 nr1int->Write(0,TObject::kOverwrite) ;
315 ar2->Write(0,TObject::kOverwrite) ;
316 br2->Write(0,TObject::kOverwrite) ;
317 cr2->Write(0,TObject::kOverwrite) ;
318 mr2->Write(0,TObject::kOverwrite) ;
319 sr2->Write(0,TObject::kOverwrite) ;
320 nr2->Write(0,TObject::kOverwrite) ;
321 nr2int->Write(0,TObject::kOverwrite) ;
322 fout.Close() ;
323
324}
325
326//-----------------------------------------------------------------------------
327Double_t PeakPosition(Double_t pt){
328 //Fit to the measured peak position
329 return 4.99292e-003*exp(-pt*9.32300e-001)+1.34944e-001 ;
330}
331//-----------------------------------------------------------------------------
332Double_t PeakWidth(Double_t pt){
333 //fit to the measured peak width
334 Double_t a=0.0068 ;
335 Double_t b=0.0025 ;
336 Double_t c=0.000319 ;
337 return TMath::Sqrt(a*a+b*b/pt/pt+c*c*pt*pt) ;
338}
339
340//-----------------------------------------------------------------------------
341Double_t CB(Double_t * x, Double_t * par){
342 //Parameterization of Real/Mixed ratio
343 Double_t m=par[1] ;
344 Double_t s=par[2] ;
345 Double_t dx=(x[0]-m)/s ;
346 return par[0]*exp(-dx*dx/2.)+par[3]+par[4]*(x[0]-kMean) ;
347}
348
349//-----------------------------------------------------------------------------
350Double_t CB2(Double_t * x, Double_t * par){
351 //Another parameterization of Real/Mixed ratio7TeV/README
352 Double_t m=par[1] ;
353 Double_t s=par[2] ;
354 Double_t dx=(x[0]-m)/s ;
355 return par[0]*exp(-dx*dx/2.)+par[3]+par[4]*(x[0]-kMean)+par[5]*(x[0]-kMean)*(x[0]-kMean) ;
356}
357//-----------------------------------------------------------------------------
358Double_t CBs(Double_t * x, Double_t * par){
359 //Parameterizatin of signal
360 Double_t m=par[1] ;
361 Double_t s=par[2] ;
362 Double_t dx=(x[0]-m)/s ;
363 return par[0]*exp(-dx*dx/2.)/TMath::Sqrt(TMath::TwoPi())/s+par[3] ;
364}
365//-----------------------------------------------------------------------------
366Double_t BG1(Double_t * x, Double_t * par){
367 //Normalizatino of Mixed
368 return par[0]+par[1]*(x[0]-kMean) ;
369}
370//-----------------------------------------------------------------------------
371Double_t BG2(Double_t * x, Double_t * par){
372 //Another normalization of Mixed
373 return par[0]+par[1]*(x[0]-kMean)+par[2]*(x[0]-kMean)*(x[0]-kMean) ;
374}
375
376
377//-----------------------------------------------------------------------------
378PPRstyle()
379{
380
381 //////////////////////////////////////////////////////////////////////
382 //
383 // ROOT style macro for the TRD TDR
384 //
385 //////////////////////////////////////////////////////////////////////
386
387 gStyle->SetPalette(1);
388 gStyle->SetCanvasBorderMode(-1);
389 gStyle->SetCanvasBorderSize(1);
390 gStyle->SetCanvasColor(10);
391
392 gStyle->SetFrameFillColor(10);
393 gStyle->SetFrameBorderSize(1);
394 gStyle->SetFrameBorderMode(-1);
395 gStyle->SetFrameLineWidth(1.2);
396 gStyle->SetFrameLineColor(1);
397
398 gStyle->SetHistFillColor(0);
399 gStyle->SetHistLineWidth(1);
400 gStyle->SetHistLineColor(1);
401
402 gStyle->SetPadColor(10);
403 gStyle->SetPadBorderSize(1);
404 gStyle->SetPadBorderMode(-1);
405
406 gStyle->SetStatColor(10);
407 gStyle->SetTitleColor(kBlack,"X");
408 gStyle->SetTitleColor(kBlack,"Y");
409
410 gStyle->SetLabelSize(0.04,"X");
411 gStyle->SetLabelSize(0.04,"Y");
412 gStyle->SetLabelSize(0.04,"Z");
413 gStyle->SetTitleSize(0.04,"X");
414 gStyle->SetTitleSize(0.04,"Y");
415 gStyle->SetTitleSize(0.04,"Z");
416 gStyle->SetTitleFont(42,"X");
417 gStyle->SetTitleFont(42,"Y");
418 gStyle->SetTitleFont(42,"X");
419 gStyle->SetLabelFont(42,"X");
420 gStyle->SetLabelFont(42,"Y");
421 gStyle->SetLabelFont(42,"Z");
422 gStyle->SetStatFont(42);
423
424 gStyle->SetTitleOffset(1.0,"X");
425 gStyle->SetTitleOffset(1.4,"Y");
426
427 gStyle->SetFillColor(kWhite);
428 gStyle->SetTitleFillColor(kWhite);
429
430 gStyle->SetOptDate(0);
431 gStyle->SetOptTitle(1);
432 gStyle->SetOptStat(0);
433 gStyle->SetOptFit(111);
434
435}
436