]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/PHOSTasks/PHOS_pp_pi0/macros/makeMmixPU_CB.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGGA / PHOSTasks / PHOS_pp_pi0 / macros / makeMmixPU_CB.C
1 void makeMmixPU_CB(const TString histoFile="LHC11a_pass4_20130913.root",
2                    const Int_t nSigma=2,
3                    const char* module="A10")
4 {
5   //Fit Real/Mixed ratio, normalize Mixed and subtract it from Real.
6   // The pi0 peak if fitted by the Crystal Ball function, 
7   // the background is fitted by pol1 or pol2
8
9   TString hMassName;
10   TFile * file = new TFile(histoFile) ;
11   THashList *hList = (THashList*)file->Get("histESD");
12   char key[125] ;
13   hMassName = "hMassPt";
14   hMassName += module;
15   TH2F * h   = (TH2F*)hList->FindObject(hMassName) ;
16
17   hMassName = "hMiMassPt";
18   hMassName += module;
19   TH2F * hm  = (TH2F*)hList->FindObject(hMassName) ;
20
21   TH1F * hev = (TH1F*)hList->FindObject("hSelEvents") ;
22
23   // Array of pt bins
24   Int_t nPadX = 6, nPadY = 3;
25   Int_t nbin=18 ;
26   Double_t xa[]={0.6,0.8,1.,1.2,1.4,1.6,1.8,2.0,2.2,2.4,2.6,3.0,3.5,4.0,5.,6.,8.,10.,12.} ;
27   PPRstyle();
28   gStyle->SetPadLeftMargin(0.16);
29   gStyle->SetPadRightMargin(0.01);
30   gStyle->SetPadTopMargin(0.02);
31   gStyle->SetPadBottomMargin(0.11);
32   gStyle->SetTitleX(0.80);
33   gStyle->SetTitleY(0.98);
34   
35   //Fit real only 
36   //Linear Bg
37   char key2[155];
38   sprintf(key,"Mix_CB") ;
39   sprintf(key2,"%s_mr1",key) ;
40   TH1D * mr1 = new TH1D(key2,"Mass",nbin,xa) ;
41   sprintf(key2,"%s_sr1",key) ;
42   TH1D * sr1 = new TH1D(key2,"Width",nbin,xa) ;
43   sprintf(key2,"%s_rms1",key) ;
44   TH1D * rms1 = new TH1D(key2,"Width",nbin,xa) ;
45   sprintf(key2,"%s_nnr1",key) ;
46   TH1D * nnr1 = new TH1D(key2,"Mass",nbin,xa) ;
47   sprintf(key2,"%s_ar1",key) ;
48   TH1D * ar1 = new TH1D(key2,"Width",nbin,xa) ;
49   sprintf(key2,"%s_yr1",key) ;
50   TH1D * nr1 = new TH1D(key2,"Raw yield",nbin,xa) ;
51   sprintf(key2,"%s_yr1int",key) ;
52   TH1D * nr1int = new TH1D(key2,"Raw yield, integrated",nbin,xa) ;
53
54   //Quadratic Bg
55   sprintf(key2,"%s_mr2",key) ;
56   TH1D * mr2 = new TH1D(key2,"Mass",nbin,xa) ;
57   sprintf(key2,"%s_sr2",key) ;
58   TH1D * sr2 = new TH1D(key2,"Width",nbin,xa) ;
59   sprintf(key2,"%s_rms2",key) ;
60   TH1D * rms2 = new TH1D(key2,"Width",nbin,xa) ;
61   sprintf(key2,"%s_nnr2",key) ;
62   TH1D * nnr2 = new TH1D(key2,"Mass",nbin,xa) ;
63   sprintf(key2,"%s_ar2",key) ;
64   TH1D * ar2 = new TH1D(key2,"Width",nbin,xa) ;
65   sprintf(key2,"%s_yr2",key) ;
66   TH1D * nr2 = new TH1D(key2,"Raw yield",nbin,xa) ;
67   sprintf(key2,"%s_yr2int",key) ;
68   TH1D * nr2int = new TH1D(key2,"Raw yield, integrated",nbin,xa) ;
69
70
71   TF1 * fit1 = new TF1("fit",CB,0.,1.,7) ;
72   fit->SetParName(0,"A") ;
73   fit->SetParName(1,"m_{0}") ;
74   fit->SetParName(2,"#sigma") ;
75   fit->SetParName(3,"a_{0}") ;
76   fit->SetParName(4,"a_{1}") ;
77   fit->SetLineWidth(2) ;
78   fit->SetLineColor(2) ;
79
80   TF1 * fgs = new TF1("gs",CBs,0.,1.,5) ;
81   fgs->SetLineColor(2) ;
82   fgs->SetLineWidth(1) ;
83
84   TF1 * fit2 = new TF1("fit2",CB2,0.,1.,8) ;
85   fit2->SetParName(0,"A") ;
86   fit2->SetParName(1,"m_{0}") ;
87   fit2->SetParName(2,"#sigma") ;
88   fit2->SetParName(3,"a_{0}") ;
89   fit2->SetParName(4,"a_{1}") ;
90   fit2->SetLineWidth(2) ;
91   fit2->SetLineColor(4) ;
92   fit2->SetLineStyle(2) ;
93
94   TF1 * fbg1 = new TF1("bg1",BG1,0.,1.,2) ;
95   TF1 * fbg2 = new TF1("bg2",BG2,0.,1.,3) ;
96
97   TCanvas * c3 = new TCanvas("mggFit1_CB_Signal","mggFitCB",10,10,1400,800) ;
98   c3->Divide(nPadX,nPadY) ;
99
100   TCanvas * cReal = new TCanvas("cReal","Mgg real events",10,10,1400,800) ;
101   cReal->Divide(nPadX,nPadY) ;
102
103   TCanvas * c1 = new TCanvas("mggFit1_CB","mggFit1",10,10,1400,800) ;
104   c1->Divide(nPadX,nPadY) ;
105   c1->cd(0) ; 
106   TCanvas * c2=0,*c4=0,*c5=0,*c6=0 ; 
107
108   TAxis * pta=h->GetYaxis() ;
109   TAxis * ma=h->GetXaxis() ;
110   for(Int_t i=1;i<=nbin;i++){
111     c1->cd(i) ;
112     Int_t imin=pta->FindBin(xa[i-1]+0.0001);
113     Int_t imax=pta->FindBin(xa[i]-0.0001) ;
114     Double_t pt=(xa[i]+xa[i-1])/2. ;
115     TH1D * hp = h->ProjectionX("re",imin,imax) ;
116     hp->Sumw2() ;
117     hp->SetNdivisions(505,"X");
118     TH1D * hpm= hm->ProjectionX("mi",imin,imax) ;
119     hpm->Sumw2() ;
120     hpm->SetNdivisions(505,"X");
121     if(i<15){
122       hp ->Rebin(2) ;
123       hpm->Rebin(2) ;
124     }
125     else{
126       hp ->Rebin(5) ;
127       hpm->Rebin(5) ;
128     }
129     hp ->SetNdivisions(506);
130     hpm->SetNdivisions(506);
131     hp ->SetLabelSize(0.05,"X");
132     hpm->SetLabelSize(0.05,"X");
133     hp ->SetTitleSize(0.00,"X");
134     hpm->SetTitleSize(0.00,"X");
135     hp ->SetLabelSize(0.05,"Y");
136     hpm->SetLabelSize(0.05,"Y");
137     hp ->SetTitleSize(0.05,"Y");
138     hpm->SetTitleSize(0.05,"Y");
139     hp ->SetTitleOffset(1.0,"X");
140     hpm->SetTitleOffset(1.0,"X");
141     if (i>12) {
142       hp ->SetLabelSize(0.05,"X");
143       hpm->SetLabelSize(0.05,"X");
144       hp ->SetTitleSize(0.05,"X");
145       hpm->SetTitleSize(0.05,"X");
146       hp ->SetXTitle("m_{#gamma#gamma} (GeV/c^{2})");
147       hpm->SetXTitle("m_{#gamma#gamma} (GeV/c^{2})");
148     }
149     // //Asign errors to the zero bins
150     // for(Int_t ib=1; ib<=hp->GetNbinsX();ib++){if(hp->GetBinContent(ib)==0)hp->SetBinError(ib,1.);}
151     // for(Int_t ib=1; ib<=hp->GetNbinsX();ib++){if(hpm->GetBinContent(ib)==0)hpm->SetBinError(ib,1.);}
152
153     TH1D * hpm2   = (TH1D*)hpm->Clone("Bg1") ;
154     TH1D * hpcopy = (TH1D*)hp ->Clone("hpcopy") ;
155     TH1D * hp2    = (TH1D*)hp ->Clone("hp2") ;
156     TH1D * hpReal = (TH1D*)hp ->Clone("hpReal") ;
157     hpcopy->Divide(hpm) ;
158     sprintf(key,"%3.1f<p_{T}<%3.1f GeV/c",xa[i-1],xa[i]) ;
159     hpcopy->SetTitle(key) ;
160     hpReal->SetTitle(key) ;
161     hpReal->SetLineWidth(2);
162     hpcopy->SetMarkerStyle(20) ;
163     hpcopy->SetMarkerSize(0.7) ;
164
165     Double_t mInit = 0.136., wInit = 0.005;
166     if(i==1) {
167       mInit = 0.139;
168       wInit = 0.007;
169     }
170     if(i==2) {
171       mInit = 0.136;
172       wInit = 0.007;
173     }
174     fit1->SetParameters(0.001,mInit,wInit,9.,0.5,0.0013,0.) ;
175     fit1->SetParLimits(2,0.003,0.010) ;
176     fit1->SetParLimits(1,0.132,0.139) ;
177     fit1->SetParLimits(0,0.,1e+6) ;
178 //    fit1->SetParLimits(3,0.,10.) ;
179 //    fit1->SetParLimits(4,0.,1e+2) ;
180     fit1->FixParameter(3,1.60) ;
181     fit1->FixParameter(4,1.27) ;
182
183     //Select fitting range
184     Double_t rangeMin=0.06 ;
185     Double_t rangeMax=0.22 ;
186     if(i>=16)rangeMax=0.3 ; //More points to fix background
187     hpcopy->Fit(fit1,"NQ","",rangeMin,rangeMax) ;
188     hpcopy->Fit(fit1,"MQ","",rangeMin,rangeMax) ;
189
190     fit2->SetParameters(fit1->GetParameters()) ;
191     fit2->SetParameter(2,wInit);
192     fit2->SetParLimits(2,0.003,0.008) ;
193     fit2->SetParLimits(1,0.130,0.142) ;
194     fit2->SetParLimits(0,0.,1e+6) ;
195 //    fit2->SetParLimits(3,0.,10.) ;
196 //    fit2->SetParLimits(4,0.,1e+2) ;
197     fit2->FixParameter(3,1.60) ;
198     fit2->FixParameter(4,1.27) ;
199     fit2->SetParameter(7,0.) ;
200
201     hpcopy->Fit(fit2,"+NQ","",rangeMin,rangeMax) ;
202     hpcopy->Fit(fit2,"+MQ","",rangeMin,rangeMax) ;
203
204     c1->cd(i) ;
205     hpcopy->SetAxisRange(0.065,0.229,"X") ;
206     hpcopy->Draw() ;
207
208     cReal->cd(i) ;
209     hpReal->SetAxisRange(0.065,0.229,"X") ;
210     hpReal->Draw("ehist") ;
211
212     fbg1->SetParameters(fit1->GetParameter(5),fit1->GetParameter(6)); 
213     fbg2->SetParameters(fit2->GetParameter(5),fit2->GetParameter(6),fit2->GetParameter(7)); 
214     Double_t intRangeMin = PeakPosition(pt)-nSigma*PeakWidth(pt) ;
215     Double_t intRangeMax = PeakPosition(pt)+nSigma*PeakWidth(pt) ;
216     Int_t    intBinMin   = hp->GetXaxis()->FindBin(intRangeMin) ;
217     Int_t    intBinMax   = hp->GetXaxis()->FindBin(intRangeMax) ;
218     Double_t errStat     = hpm->Integral(intBinMin,intBinMax); 
219
220     hpm ->Multiply(fbg1) ;
221     hpm2->Multiply(fbg2) ;
222     hp  ->Add(hpm,-1.) ;
223     hp2 ->Add(hpm2,-1.) ;
224
225     c3->cd(i) ;
226
227     if(i<15)
228       fgs->SetParameters(hp->Integral(32,36)/5.,fit1->GetParameter(1),fit1->GetParameter(2),fit1->GetParameter(3),fit1->GetParameter(4)) ;
229     else
230       fgs->SetParameters(hp->Integral(13,15)/3.,fit1->GetParameter(1),fit1->GetParameter(2),fit1->GetParameter(3),fit1->GetParameter(4)) ;
231 //      fgs->SetParameters(hp->Integral(13,15)/3.,0.135,0.008,0.7,2.) ;
232
233     fgs->SetParLimits(0,0.,1e+06) ;
234     fgs->SetParLimits(1,0.130,0.142) ;
235     fgs->SetParLimits(2,0.003,0.008) ;
236     fgs->FixParameter(3,1.60) ;
237     fgs->FixParameter(4,1.27) ;
238     
239     hp->Fit(fgs,"Q","",rangeMin,rangeMax) ;   
240     hp->SetMaximum(hp2->GetMaximum()*1.1) ;
241     hp->SetMinimum(hp2->GetMinimum()*1.1) ;
242     hp->SetMarkerStyle(20) ;
243     hp->SetMarkerSize(0.7) ;
244     mr1->SetBinContent (i,fgs->GetParameter(1)) ;
245     mr1->SetBinError   (i,fgs->GetParError(1)) ;
246     sr1->SetBinContent (i,TMath::Abs(fgs->GetParameter(2))) ;
247     sr1->SetBinError   (i,fgs->GetParError(2)) ;
248     nnr1->SetBinContent(i,fgs->GetParameter(3)) ;
249     nnr1->SetBinError  (i,fgs->GetParError(3)) ;
250     ar1->SetBinContent (i,fgs->GetParameter(4)) ;
251     ar1->SetBinError   (i,fgs->GetParError(4)) ;
252
253     Double_t y=fgs->Integral(intRangeMin,intRangeMax)/hp->GetXaxis()->GetBinWidth(1) ;
254     nr1->SetBinContent(i,y) ;
255     if(y>0)
256       rms1->SetBinContent(i,fgs->CentralMoment(2.,0.05,0.2)) ;
257     Double_t ey=0 ;
258     if(fgs->GetParameter(0)!=0. && fgs->GetParameter(2)!=0.){
259       Double_t en=fgs->GetParError(0)/fgs->GetParameter(0) ;
260       Double_t es=fgs->GetParError(2)/fgs->GetParameter(2) ;
261       ey=y*TMath::Sqrt(en*en+es*es) ;
262     }
263     nr1->SetBinError(i,ey) ;
264
265     Double_t npiInt=hp->Integral(intBinMin,intBinMax) ;
266     Double_t norm=fbg1->GetParameter(0) ;
267     Double_t normErr=fbg1->GetParError(0) ;
268     if(npiInt>0.){
269       nr1int->SetBinContent(i,npiInt) ;
270       nr1int->SetBinError(i,TMath::Sqrt(npiInt + norm*errStat + normErr*normErr*errStat*errStat + norm*norm*errStat)) ;
271     }
272
273     hp2->GetXaxis()->SetRangeUser(0.05,0.25) ;
274     hp2->SetMaximum(hp2->GetMaximum()*1.1) ;
275     hp2->SetMinimum(hp2->GetMinimum()*1.1) ;
276     hp2->SetMarkerStyle(20) ;
277     hp2->SetMarkerSize(0.7) ;
278
279     hp2->Fit(fgs,"Q","",rangeMin,rangeMax) ;
280     mr2->SetBinContent (i,fgs->GetParameter(1)) ;
281     mr2->SetBinError   (i,fgs->GetParError(1)) ;
282     sr2->SetBinContent (i,TMath::Abs(fgs->GetParameter(2))) ;
283     sr2->SetBinError   (i,fgs->GetParError(2)) ;
284     nnr2->SetBinContent(i,fgs->GetParameter(3)) ;
285     nnr2->SetBinError  (i,fgs->GetParError(3)) ;
286     ar2->SetBinContent (i,fgs->GetParameter(4)) ;
287     ar2->SetBinError   (i,fgs->GetParError(4)) ;
288
289     y=fgs->Integral(intRangeMin,intRangeMax)/hp->GetXaxis()->GetBinWidth(1) ;
290     nr2->SetBinContent(i,y) ;
291     ey=0 ;
292     if(fgs->GetParameter(0)!=0. && fgs->GetParameter(2)!=0.){
293       Double_t en=fgs->GetParError(0)/fgs->GetParameter(0) ;
294       Double_t es=fgs->GetParError(2)/fgs->GetParameter(2) ;
295       ey=y*TMath::Sqrt(en*en+es*es) ;
296     }
297     nr2->SetBinError(i,ey) ;
298     npiInt=hp2->Integral(intBinMin,intBinMax) ;
299     norm=fbg2->GetParameter(0) ;
300     normErr=fbg2->GetParError(0) ;
301     if(npiInt>0.){
302       nr2int->SetBinContent(i,npiInt) ;
303       nr2int->SetBinError(i,TMath::Sqrt(npiInt + norm*errStat + normErr*normErr*errStat*errStat + norm*norm*errStat)) ;
304     } 
305
306     hp2->SetTitle(key) ;
307     hp2->SetAxisRange(0.065,0.229,"X") ;
308     hp2->Draw() ;
309     delete hp ;
310 //    delete hp2 ;
311 //    delete hpcopy ;
312     delete hpm ;
313     delete hpm2 ;
314     c1->Update() ;
315     c3->Update() ;
316     cReal->Update() ;
317   }
318
319   if (c3) c3->Print("Pi0_Signal_CB.eps") ;
320   if (c1) c1->Print("Pi0_Ratio_CB.eps") ;
321   if (cReal) cReal->Print("Pi0_InvMass_CB.eps") ;
322
323   //Normalize by the number of non-pileup events
324   Double_t nMBOR   = hev->GetBinContent(2); // MBOR events pass4
325   // Double_t nPU     = hev->GetBinContent(7);
326   // Double_t nMBOR   = hev->GetBinContent(1); // MBOR events pass3
327   // Double_t nPU     = hev->GetBinContent(6);
328   Double_t nevents = nMBOR;
329   // Double_t nevents = nMBOR - nPU;
330   printf("==============\nN events = %d\n==============\n",nevents);
331
332   nr1   ->Scale(1./nevents) ;
333   nr1int->Scale(1./nevents) ;
334   nr2   ->Scale(1./nevents) ;
335   nr2int->Scale(1./nevents) ;
336
337
338   TFile fout("LHC11a_FitResult.root","update");
339   mr1   ->Write(0,TObject::kOverwrite) ;
340   sr1   ->Write(0,TObject::kOverwrite) ;
341   rms1  ->Write(0,TObject::kOverwrite) ;
342   nnr1  ->Write(0,TObject::kOverwrite) ;
343   ar1   ->Write(0,TObject::kOverwrite) ;
344   nr1   ->Write(0,TObject::kOverwrite) ;
345   nr1int->Write(0,TObject::kOverwrite) ;
346   mr2   ->Write(0,TObject::kOverwrite) ;
347   sr2   ->Write(0,TObject::kOverwrite) ;
348   rms2  ->Write(0,TObject::kOverwrite) ;
349   nnr2  ->Write(0,TObject::kOverwrite) ;
350   ar2   ->Write(0,TObject::kOverwrite) ;
351   nr2   ->Write(0,TObject::kOverwrite) ;
352   nr2int->Write(0,TObject::kOverwrite) ;
353   fout.Close() ;
354
355 }
356
357
358 //-----------------------------------------------------------------------------
359 const Double_t kMean=0.135 ; //Approximate peak position to facilitate error estimate
360
361 Double_t PeakPosition(Double_t pt){
362   //Fit to the measured peak position
363   // return 4.99292e-003*exp(-pt*9.32300e-001)+1.34944e-001 ;
364   return 0.57*exp(-pt*7.62)+0.13592 ;
365 }
366 Double_t PeakWidth(Double_t pt){
367   //fit to the measured peak width
368   return 1.60935e-02*exp(-pt*2.25609e+00)+4.65743e-03 ;
369   // return 0.0068 ;
370 }
371  
372 Double_t CB(Double_t * x, Double_t * par){
373   Double_t m=par[1] ;
374   Double_t s=par[2] ;
375   Double_t n=par[3] ;
376   Double_t a=par[4] ;
377   Double_t dx=(x[0]-m)/s ;
378   if(dx>-a)
379     return par[0]*exp(-dx*dx/2.)+par[5]+par[6]*(x[0]-kMean) ;
380   else{
381     Double_t A=TMath::Power((n/TMath::Abs(a)),n)*TMath::Exp(-a*a/2) ;
382     Double_t B=n/TMath::Abs(a)-TMath::Abs(a) ;
383     return par[0]*A*TMath::Power((B-dx),-n)+par[5]+par[6]*(x[0]-kMean) ;
384   }
385
386 }
387 Double_t CB2(Double_t * x, Double_t * par){
388   Double_t m=par[1] ;
389   Double_t s=par[2] ;
390   Double_t n=par[3] ;
391   Double_t a=par[4] ;
392   Double_t dx=(x[0]-m)/s ;
393   if(dx>-a)
394     return par[0]*exp(-dx*dx/2.)+par[5]+par[6]*(x[0]-kMean)+par[7]*(x[0]-kMean)*(x[0]-kMean) ;
395   else{
396     Double_t A=TMath::Power((n/TMath::Abs(a)),n)*TMath::Exp(-a*a/2) ;
397     Double_t B=n/TMath::Abs(a)-TMath::Abs(a) ;
398     return par[0]*A*TMath::Power((B-dx),-n)+par[5]+par[6]*(x[0]-kMean)+par[7]*(x[0]-kMean)*(x[0]-kMean) ;
399   }
400
401 }
402 Double_t CBs(Double_t * x, Double_t * par){
403   Double_t m=par[1] ;
404   Double_t s=par[2] ;
405   Double_t n=par[3] ;
406   Double_t a=par[4] ;
407   Double_t dx=(x[0]-m)/s ;
408   if(dx>-a)
409     return par[0]*exp(-dx*dx/2.) ;
410   else{
411     Double_t A=TMath::Power((n/TMath::Abs(a)),n)*TMath::Exp(-a*a/2) ;
412     Double_t B=n/TMath::Abs(a)-TMath::Abs(a) ;
413     return par[0]*A*TMath::Power((B-dx),-n) ;
414   }
415 }
416 Double_t BG1(Double_t * x, Double_t * par){
417   //Normalizatino of Mixed
418   return par[0]+par[1]*(x[0]-kMean) ;
419 }
420 Double_t BG2(Double_t * x, Double_t * par){
421   //Another normalization of  Mixed
422   return par[0]+par[1]*(x[0]-kMean)+par[2]*(x[0]-kMean)*(x[0]-kMean) ;
423 }
424
425 //-----------------------------------------------------------------------------
426 PPRstyle()
427 {
428
429   //////////////////////////////////////////////////////////////////////
430   //
431   // ROOT style macro for the TRD TDR
432   //
433   //////////////////////////////////////////////////////////////////////
434
435   gStyle->SetPalette(1);
436   gStyle->SetCanvasBorderMode(-1);
437   gStyle->SetCanvasBorderSize(1);
438   gStyle->SetCanvasColor(10);
439
440   gStyle->SetFrameFillColor(10);
441   gStyle->SetFrameBorderSize(1);
442   gStyle->SetFrameBorderMode(-1);
443   gStyle->SetFrameLineWidth(1.2);
444   gStyle->SetFrameLineColor(1);
445
446   gStyle->SetHistFillColor(0);
447   gStyle->SetHistLineWidth(1);
448   gStyle->SetHistLineColor(1);
449
450   gStyle->SetPadColor(10);
451   gStyle->SetPadBorderSize(1);
452   gStyle->SetPadBorderMode(-1);
453
454   gStyle->SetStatColor(10);
455   gStyle->SetTitleColor(kBlack,"X");
456   gStyle->SetTitleColor(kBlack,"Y");
457
458   gStyle->SetLabelSize(0.04,"X");
459   gStyle->SetLabelSize(0.04,"Y");
460   gStyle->SetLabelSize(0.04,"Z");
461   gStyle->SetTitleSize(0.04,"X");
462   gStyle->SetTitleSize(0.04,"Y");
463   gStyle->SetTitleSize(0.04,"Z");
464   gStyle->SetTitleFont(42,"X");
465   gStyle->SetTitleFont(42,"Y");
466   gStyle->SetTitleFont(42,"X");
467   gStyle->SetLabelFont(42,"X");
468   gStyle->SetLabelFont(42,"Y");
469   gStyle->SetLabelFont(42,"Z");
470   gStyle->SetStatFont(42);
471
472   gStyle->SetTitleOffset(1.0,"X");
473   gStyle->SetTitleOffset(1.4,"Y");
474
475   gStyle->SetFillColor(kWhite);
476   gStyle->SetTitleFillColor(kWhite);
477
478   gStyle->SetOptDate(0);
479   gStyle->SetOptTitle(1);
480   gStyle->SetOptStat(0);
481   gStyle->SetOptFit(0);
482
483 }
484