]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/EVCHAR/correctSPDdNdEtapp.C
Add unicor Kink evchar
[u/mrichter/AliRoot.git] / PWG2 / EVCHAR / correctSPDdNdEtapp.C
CommitLineData
bb5d72d4 1/*************************************************************************
2* Macro correctSPDdNdEtapp *
3* To read data and correction histograms produced with *
4* AliAnalysisTaskSPDdNdEta and apply corrections to data *
5* *
6* Author: M. Nicassio (INFN Bari) *
7* Contact: Maria.Nicassio@ba.infn.it, Domenico.Elia@ba.infn.it *
8**************************************************************************/
9
10#include "Riostream.h"
11#include "TROOT.h"
12#include "TStyle.h"
13#include "TMath.h"
14#include "TH1.h"
15#include "TH2.h"
16#include "TFile.h"
17#include "TList.h"
18#include "TTree.h"
19#include "TLegend.h"
20#include "TLegendEntry.h"
21#include "TCanvas.h"
22
23void correctSPDdNdEtapp (Char_t* fileRaw, Char_t* fileCorr, Char_t* fileMC, Char_t* fileAccPb,
24 Bool_t kPrimariesTest, Bool_t kAccPb,
25 Int_t binVtxStart, Int_t binVtxStop,
26 Double_t cutAccPb, Double_t cutAcc, Double_t cutBkg,
27 Double_t cutAlgEff, Double_t cutDec, Double_t cutVtx) {
28
29//lim vtx: 5-15 ->[-10cm,10cm[; 0-20 ->[-20cm,20cm[
30
31 gStyle->SetOptLogy(kFALSE);
32 gStyle->SetFrameLineWidth(2);
33 gStyle->SetCanvasColor(0);
34 gStyle->SetPadColor(0);
35 gStyle->SetHistLineWidth(2);
36 gStyle->SetLabelSize(0.05, "x");
37 gStyle->SetLabelSize(0.05, "y");
38 gStyle->SetTitleSize(0.05, "x");
39 gStyle->SetTitleSize(0.05, "y");
40 gStyle->SetTitleOffset(1.1, "x");
41 gStyle->SetPadBottomMargin(0.14);
42// gStyle->SetPalette(1,0);
43 gStyle->SetFillColor(0);
44
45
46 //Input files
47 TFile *f_accPb = new TFile(fileAccPb);
48 TFile *f_corr = new TFile(fileCorr);
49 TFile *f_MC = new TFile(fileMC);
50 TFile *f_raw = new TFile(fileRaw);
51
52 TList *l_corr = (TList*)f_corr->Get("cOutput");
53 TList *l_MC = (TList*)f_MC->Get("cOutput");
54 TList *l_raw = (TList*)f_raw->Get("cOutput");
55
56 //Corrections
57 //Background
58 TH2F *hBackgroundCorr = new TH2F("backgroundCorr","Background correction",60,-3.,3.,20,-20.,20.);
59 hBackgroundCorr->Sumw2();
60 TH2F *hBackgroundCorrNum = (TH2F*)l_corr->At(31);
61
62 Int_t corrHistoName;
63 if (kPrimariesTest) {
64 corrHistoName = 31;
65 } else {
66 corrHistoName = 32;
67 }
68
69 TH2F *hBackgroundCorrDen = (TH2F*)l_corr->At(corrHistoName);
70 hBackgroundCorr->Divide(hBackgroundCorrNum,hBackgroundCorrDen,1.,1.);
71
72 //Algorithm efficiency
73 TH2F *hTrackletRecEff = new TH2F("trackletRecEffCorr","Tracklet algorithm efficiency",60,-3,3,20,-20,20);
74 hTrackletRecEff->Sumw2();
75 TH2F *hTrackletEffCorr_num = (TH2F*)l_corr->At(33);
76 hTrackletRecEff->Divide(hBackgroundCorrNum,hTrackletEffCorr_num,1.,1.); // 1/Efficiency=correction_weight
77
78 //Correction for particles that do not reach the sensitive layer
79 TH2F *hDecPartLay2Corr = new TH2F("decPartLay2Corr","",60,-3.,3.,20,-20.,20.);
80 hDecPartLay2Corr->Sumw2();
81 TH2F *hPrimaries_evtTrVtx = (TH2F*) l_corr->At(34);
82 TH2F *hDetectablePrimariesLay2_evtTrVtx = (TH2F*) l_corr->At(35);
83 hDecPartLay2Corr->Divide(hPrimaries_evtTrVtx,hDetectablePrimariesLay2_evtTrVtx,1.,1.);
84
85 //Acceptance
86 TH2F *hAccCorr = new TH2F("accEffCorr_test","",60,-3,3,20,-20,20);
87 hAccCorr->Sumw2();
88 hAccCorr->Divide(hTrackletEffCorr_num,hDetectablePrimariesLay2_evtTrVtx,1.,1.);
89
90 TH2F *hAccCorrPb = (TH2F*) f_accPb->Get("AccTrac");
91
92 //Vertex-trigger correction
93 TH2F *hVertexTrigCorrEtaNum = (TH2F*)l_corr->At(36);
94 TH2F *hTriggerCorrEta_norm = (TH2F*)l_corr->At(37);
95 TH2F *hVertexCorrEta_norm = (TH2F*)l_corr->At(34);
96
97 TH2F *hVertexCorrEta = new TH2F("vertexCorrEta","",60,-3.,3.,20,-20.,20.);
98 hVertexCorrEta->Sumw2();
99 hVertexCorrEta->Divide(hTriggerCorrEta_norm,hVertexCorrEta_norm,1.,1.);
100
101 TH2F *hTriggerCorrEta = new TH2F("triggerEtaCorr","",60,-3.,3.,20,-20.,20.);
102 hTriggerCorrEta->Sumw2();
103 hTriggerCorrEta->Divide(hVertexTrigCorrEtaNum,hTriggerCorrEta_norm,1.,1.);
104
105 TH2F *hVertexTrigCorrEta = new TH2F("vertexTrigEtaCorr","",60,-3.,3.,20,-20.,20.);
106 hVertexTrigCorrEta->Sumw2();
107 hVertexTrigCorrEta->Divide(hVertexTrigCorrEtaNum,hVertexCorrEta_norm,1.,1.);
108
109 TH2F *hVertexTrigCorrNum = (TH2F*)l_corr->At(38);
110 TH2F *hTriggerCorr_norm = (TH2F*)l_corr->At(40);
111 TH2F *hVertexCorr_norm = (TH2F*)l_corr->At(39);
112
113 TH2F *hVertexCorr = new TH2F("vertexCorr","",200,0.,200.,20,-20.,20.);
114 hVertexCorr->Sumw2();
115 hVertexCorr->Divide(hTriggerCorr_norm,hVertexCorr_norm,1.,1.);
116
117 TH2F *hTriggerCorr = new TH2F("triggerCorr","",200,0.,200.,20,-20.,20.);
118 hTriggerCorr->Sumw2();
119 hTriggerCorr->Divide(hVertexTrigCorrNum,hTriggerCorr_norm,1.,1.);
120
121 TH2F *hVertexTrigCorr = new TH2F("vertexTrigCorr","",200,0.,200.,20,-20.,20.);
122 hVertexTrigCorr->Sumw2();
123 hVertexTrigCorr->Divide(hVertexTrigCorrNum,hVertexCorr_norm,1.,1.);
124
125 TH2F *hVertexMC2D = (TH2F*)l_MC->At(38);
126 TH1D *hMCvertex = new TH1D(*(hVertexMC2D->ProjectionY("MCvertex",0,-1,"e"))); //MC vertex distrib all events
127// hMCvertex->Draw();
128
129 //Raw distributions
130
131 //...track level
132 TH2F *hSPDEta_2Draw_accBin = new TH2F(*((TH2F*)l_raw->At(2)));
133 hSPDEta_2Draw_accBin->SetNameTitle("SPDEta_2Draw_accBin","");
134 TH2F *hSPDEta_2Draw = new TH2F(*((TH2F*)l_raw->At(2)));
135 hSPDEta_2Draw->SetNameTitle("SPDEta_2Draw","");
136
137 //...event level
138 TH2F *hESDVertexvsNtracklets = (TH2F*)l_raw->At(0);
139 TH2F *hESDVertexvsNtracklets_trigEvts = (TH2F*)l_raw->At(1);
140
141 TH2F *hESDVertexvsNtracklets_Corr = new TH2F("ESDVertexvsNtracklets_Corr","",200,0.,200.,20,-20.,20.);
142
143 TH1D *hSPDvertex = (TH1D*)l_raw->At(18);
144 if (kPrimariesTest) {
145 hSPDEta_2Draw = new TH2F(*((TH2F*)l_MC->At(31)));
146 hSPDEta_2Draw->SetNameTitle("SPDEta_2Draw","");
147 hESDVertexvsNtracklets = new TH2F(*((TH2F*)l_MC->At(39)));
148 hESDVertexvsNtracklets_trigEvts = new TH2F(*((TH2F*)l_MC->At(40)));
149
150 hSPDvertex = hVertexCorr_norm->ProjectionY("MCvertex_RV",0,-1,"e");
151 }
152
153 if (!kPrimariesTest) hSPDEta_2Draw->Rebin2D(2,2,"");
154
155 //Vertex correction
156
157 //Computing the percentage of triggered events with vertex and tracklet mult!=0 for each rec vertex bin
158 TH1F* fHistZdistribTrigVtxEvt = new TH1F("fHistZdistribTrigVtxEvt", "",20,-20.,20.);
159 fHistZdistribTrigVtxEvt->Sumw2();
160 for (Int_t ivtx=binVtxStart; ivtx<binVtxStop; ivtx++) {
161
162 Float_t fracZ = (hESDVertexvsNtracklets->ProjectionY("_y",2,200,"e")->GetBinContent(ivtx+1))/Float_t (hESDVertexvsNtracklets->ProjectionY("_y",2,200,"e")->Integral(binVtxStart+1,binVtxStop));
163 fHistZdistribTrigVtxEvt->SetBinContent(ivtx+1,fracZ);
164 }
165
166 //Computing the factor to correct the previous %
167 TH1F* fHistEvtCorrFactor = new TH1F("fHistEvtCorrFactor", "",20,-20.,20.);
168 fHistEvtCorrFactor->Sumw2();
169 for (Int_t ivtx=0; ivtx<20; ivtx++) {
170 Float_t fracZnum = (hTriggerCorr_norm->GetBinContent(1,ivtx+1))/Float_t(hTriggerCorr_norm->ProjectionY("_y",1,1,"e")->Integral());
171 Float_t fracZden = (hVertexCorr_norm->ProjectionY("_y",2,200,"e")->GetBinContent(ivtx+1))/Float_t(hVertexCorr_norm->ProjectionY("_y",2,200,"e")->Integral());
172 Float_t ratio =0.;
173 if (fracZden!=0.) ratio = fracZnum/fracZden;
174 if(ratio!=0.) fHistEvtCorrFactor->SetBinContent(ivtx+1,ratio);
175 }
176
177 TH1F* hprod = new TH1F("prod", "",20,-20.,20.);
178 hprod->Sumw2();
179
180 hprod->Multiply(fHistEvtCorrFactor,fHistZdistribTrigVtxEvt,1.,1.);
181
182 //MC distributions
183
184 Int_t histoNameRV;
185 if (kPrimariesTest) {
186 histoNameRV = 41;
187 } else {
188 histoNameRV = 42;
189 }
190
191 TH2F *Eta_RV = (TH2F*)l_MC->At(histoNameRV);
192 TH2F *Eta = (TH2F*)l_MC->At(43);
193// Eta->Draw();
194
195 TH1D *hdNdEta_RV = new TH1D("dNdEta_RV","Pseudorapidity ",60,-3,3);
196 hdNdEta_RV=Eta_RV->ProjectionX("dNdEta_RV",binVtxStart+1,binVtxStop,"e");
197 Double_t nEvtsRec=0;
198 for (Int_t ivtx=binVtxStart; ivtx<binVtxStop; ivtx++) {
199 nEvtsRec += hSPDvertex->GetBinContent(ivtx+1);
200 }
201 hdNdEta_RV->Scale(10./nEvtsRec);
202
203 TH1D *hdNdEta_Trigg = new TH1D("dNdEta_Trigg","Pseudorapidity ",60,-3,3);
204 TH2F *hEta_Trigg = new TH2F(*((TH2F*)l_corr->At(37)));
205 hdNdEta_Trigg=hEta_Trigg->ProjectionX("dNdEta_Trigg",0,-1,"e");
206 Double_t nEvtsTrigg=0;
207 nEvtsTrigg =hTriggerCorr_norm->Integral();
208 hdNdEta_Trigg->Scale(10./nEvtsTrigg);
209
210 TH1D *hdNdEta = new TH1D("dNdEta","Pseudorapidity ",60,-3.,3.);
211 hdNdEta->Sumw2();
212 hdNdEta = Eta->ProjectionX("dNdEta",0,-1,"e");
213 Double_t nEvtsTot=0;
214 for (Int_t ivtx=0; ivtx<20; ivtx++) {
215 nEvtsTot += hMCvertex->GetBinContent(ivtx+1);
216 }
217 cout<<"Number of generated events: "<<nEvtsTot <<endl;
218
219 hdNdEta->Scale(10./nEvtsTot);
220 hdNdEta->GetXaxis()->SetTitle("Pseudorapidity #eta");
221 hdNdEta->GetYaxis()->SetTitle("dN_{ch}/d#eta");
222// hdNdEta->Draw();
223
224 //Corrected distributions
225 TH2F *hSPDEta_2D_bkgCorr = new TH2F("SPDEta_2D_bkgCorr", "",60,-3.,3.,20,-20.,20.);
226 hSPDEta_2D_bkgCorr->Sumw2();
227 TH2F *hSPDEta_2D_bkgEffCorr = new TH2F("SPDEta_2D_bkgEffCorr", "",60,-3.,3.,20,-20.,20.);
228 hSPDEta_2D_bkgEffCorr->Sumw2();
229 TH2F *hSPDEta_2D_bkgEffAccCorr = new TH2F("SPDEta_2D_bkgEffAccCorr", "",120,-3.,3.,40,-20.,20.);
230 hSPDEta_2D_bkgEffAccCorr->Sumw2();
231 TH2F *hSPDEta_2D_fullyCorr = new TH2F("SPDEta_2D_fullyCorr", "",60,-3.,3.,20,-20.,20.);
232 hSPDEta_2D_fullyCorr->Sumw2();
233 TH2F *hSPDEta_2D_triggeredEvts = new TH2F("SPDEta_2D_triggeredEvts","",60,-3.,3.,20,-20.,20.);
234 hSPDEta_2D_triggeredEvts->Sumw2();
235 TH2F *hSPDEta_2D_allEvts = new TH2F("SPDEta_2D_allEvts","",60,-3.,3.,20,-20.,20.);
236 hSPDEta_2D_allEvts->Sumw2();
237
238 TH1F *hnorm_bkgCorr = new TH1F("norm_bkgCorr", "",60,-3.,3.);
239 hnorm_bkgCorr->Sumw2();
240 TH1F *hnorm_bkgEffCorr = new TH1F("norm_bkgEffCorr", "",60,-3.,3.);
241 hnorm_bkgEffCorr->Sumw2();
242 TH1F *hnorm_bkgEffAccCorr = new TH1F("norm_bkgEffAccCorr", "",60,-3.,3.);
243 hnorm_bkgEffAccCorr->Sumw2();
244 TH1F *hnorm_fullyCorr = new TH1F("norm_fullyCorr", "",60,-3.,3.);
245 hnorm_fullyCorr->Sumw2();
246 TH1F *hnorm_triggered = new TH1F("norm_triggered", "",60,-3.,3.);
247 hnorm_triggered->Sumw2();
248 TH1F *hnorm = new TH1F("norm", "",60,-3.,3.);
249 hnorm->Sumw2();
250 TH1F *hnorm_gen = new TH1F("norm_gen", "",60,-3.,3.);
251 hnorm_gen->Sumw2();
252
253 //dN/deta
254 TH1F *hdNdEta_raw = new TH1F("dNdEta_raw","", 60,-3.,3.);
255 hdNdEta_raw->Sumw2();
256 hdNdEta_raw->Add(hSPDEta_2Draw->ProjectionX("SPDEta_2Draw_x",binVtxStart+1,binVtxStop,"e"));
257 hdNdEta_raw->Scale(10./nEvtsRec);
258
259 TH1F *hdNdEta_bkgCorr = new TH1F("dNdEta_bkgCorr","", 60,-3.,3.);
260 hdNdEta_bkgCorr->Sumw2();
261 TH1F *hdNdEta_bkgEffCorr = new TH1F("dNdEta_bkgEffCorr","", 60,-3.,3.);
262 hdNdEta_bkgEffCorr->Sumw2();
263 TH1F *hdNdEta_bkgEffAccCorr = new TH1F("dNdEta_bkgEffAccCorr","", 60,-3.,3.);
264 hdNdEta_bkgEffAccCorr->Sumw2();
265 TH1F *hdNdEta_fullyCorr = new TH1F("dNdEta_fullyCorr","", 60,-3.,3.);
266 hdNdEta_fullyCorr->Sumw2();
267 TH1F *hdNdEta_triggeredEvts = new TH1F("dNdEta_triggeredEvts","", 60,-3.,3.);
268 hdNdEta_triggeredEvts->Sumw2();
269 TH1F *hdNdEta_allEvts = new TH1F("dNdEta_allEvts","", 60,-3.,3.);
270 hdNdEta_allEvts->Sumw2();
271
272 //Cutting either on weights or on relative errors of weights
273 // Acceptance Pb-Pb
274 TH2F *hAccErrorsvsAccTr = new TH2F("AccErrorsvsAccTr","",20000,0.,10000.,50,0.,1.);
275 for (Int_t jeta=0; jeta<120; jeta++) {
276 for (Int_t kvtx=0; kvtx<40; kvtx++) {
277 if (hAccCorrPb->GetBinContent(jeta+1,kvtx+1) < cutAccPb) {
278 hAccCorrPb->SetBinContent(jeta+1,kvtx+1,0.);
279 hAccCorrPb->SetBinError(jeta+1,kvtx+1,0.);
280 }
281 }
282 }
283
284 for (Int_t jeta=0; jeta<60; jeta++) {
285 for (Int_t kvtx=0; kvtx<20; kvtx++) {
286 // Background
287 Float_t relBinErr = hBackgroundCorr->GetBinError(jeta+1,kvtx+1)/hBackgroundCorr->GetBinContent(jeta+1,kvtx+1);
288 if (relBinErr>cutBkg) {
289 hBackgroundCorr->SetBinContent(jeta+1,kvtx+1,0.);
290 hBackgroundCorr->SetBinError(jeta+1,kvtx+1,0.);
291 }
292 //Alg eff
293 relBinErr = hTrackletRecEff->GetBinError(jeta+1,kvtx+1)/hTrackletRecEff->GetBinContent(jeta+1,kvtx+1);
294 if (relBinErr>cutAlgEff) {
295 hTrackletRecEff->SetBinContent(jeta+1,kvtx+1,0.);
296 hTrackletRecEff->SetBinError(jeta+1,kvtx+1,0.);
297 }
298 //Acceptance p-p
299 relBinErr = hAccCorr->GetBinError(jeta+1,kvtx+1)/hAccCorr->GetBinContent(jeta+1,kvtx+1);
300 if (relBinErr>cutAcc) {
301// if (hAccCorr->GetBinContent(jeta+1,kvtx+1)<cutAcc) {
302 hAccCorr->SetBinContent(jeta+1,kvtx+1,0.);
303 hAccCorr->SetBinError(jeta+1,kvtx+1,0.);
304 }
305 //Dec particles
306 relBinErr = hDecPartLay2Corr->GetBinError(jeta+1,kvtx+1)/hDecPartLay2Corr->GetBinContent(jeta+1,kvtx+1);
307 if (relBinErr>cutDec) {
308 hDecPartLay2Corr->SetBinContent(jeta+1,kvtx+1,0.);
309 hDecPartLay2Corr->SetBinError(jeta+1,kvtx+1,0.);
310 }
311 //Vtx Corr
312 relBinErr = hVertexCorrEta->GetBinError(jeta+1,kvtx+1)/hVertexCorrEta->GetBinContent(jeta+1,kvtx+1);
313 if (relBinErr>cutVtx) {
314 hVertexCorrEta->SetBinContent(jeta+1,kvtx+1,0.);
315 hVertexCorrEta->SetBinError(jeta+1,kvtx+1,0.);
316 }
317 }
318 }
319
320 for (Int_t ieta=0; ieta<60; ieta++) {
321 for (Int_t izeta=0; izeta<20; izeta++) {
322 if (hAccCorr->GetBinContent(ieta+1,izeta+1)==0)
323 hDecPartLay2Corr->SetBinContent(ieta+1,izeta+1,0.);
324 }
325 }
326
327 for (Int_t ieta=0; ieta<60; ieta++) {
328 for (Int_t izeta=0; izeta<20; izeta++) {
329 if (hAccCorr->GetBinContent(ieta+1,izeta+1)==0)
330 hVertexTrigCorrEta->SetBinContent(ieta+1,izeta+1,0.);
331 }
332 }
333
334 //Apply corrections at
335 //...track level
336 hSPDEta_2D_bkgCorr->Add(hSPDEta_2Draw);
337 hSPDEta_2D_bkgCorr->Multiply(hBackgroundCorr);
338 hSPDEta_2D_bkgEffCorr->Divide(hSPDEta_2D_bkgCorr,hTrackletRecEff,1.,1.);
339
340 if (kAccPb) {
341 hSPDEta_2D_bkgEffAccCorr->Divide(hSPDEta_2Draw_accBin,hAccCorrPb,1.,1.);
342 hSPDEta_2D_bkgEffAccCorr->Rebin2D(2,2,"");
343 } else {
344 hSPDEta_2D_bkgEffAccCorr->Rebin2D(2,2,"");
345 hSPDEta_2D_bkgEffAccCorr->Add(hSPDEta_2Draw);
346 hSPDEta_2D_bkgEffAccCorr->Divide(hAccCorr);
347 }
348 hSPDEta_2D_bkgEffAccCorr->Multiply(hBackgroundCorr);
349 hSPDEta_2D_bkgEffAccCorr->Divide(hTrackletRecEff);
350 hSPDEta_2D_fullyCorr->Add(hSPDEta_2D_bkgEffAccCorr);
351 hSPDEta_2D_fullyCorr->Multiply(hDecPartLay2Corr);
352 hSPDEta_2D_triggeredEvts->Add(hSPDEta_2D_fullyCorr);
353 hSPDEta_2D_triggeredEvts->Multiply(hVertexCorrEta);
354 hSPDEta_2D_allEvts->Add(hSPDEta_2D_triggeredEvts);
355 hSPDEta_2D_allEvts->Multiply(hTriggerCorrEta);
356
357 //...event level
358 Double_t evtTrigNoTracklets = hESDVertexvsNtracklets_trigEvts->ProjectionY("_y",1,1,"e")->Integral();
359 for (Int_t ivtx=0;ivtx<20;ivtx++) {
360 hESDVertexvsNtracklets_trigEvts->SetBinContent(1,ivtx+1,(hprod->GetBinContent(ivtx+1))*(evtTrigNoTracklets));
361 }
362
363 hESDVertexvsNtracklets_Corr->Add(hESDVertexvsNtracklets_trigEvts);
364 hESDVertexvsNtracklets_Corr->Multiply(hTriggerCorr);
365
366 //Filling normalization histograms
367
368 TH1D *hSPDvertexCorr_y = new TH1D("SPDvertexCorr_y","",20,-20,20);
369 hSPDvertexCorr_y = hESDVertexvsNtracklets_Corr->ProjectionY("SPDvertexCorr_y",0,-1,"e");
370 cout<<"Number of events (tracklets corrected)"<<hSPDvertexCorr_y->Integral(binVtxStart+1,binVtxStop)<<endl;
371 for (Int_t ivtx=binVtxStart; ivtx<binVtxStop; ivtx++) {
372 Double_t nEvtsivtx = hSPDvertex->GetBinContent(ivtx+1);
373 Double_t nEvtsivtxCorr = hSPDvertexCorr_y->GetBinContent(ivtx+1);
374 Double_t nEvtsivtxCorr_triggered = hESDVertexvsNtracklets_trigEvts->ProjectionY("_y",0,-1,"e")->GetBinContent(ivtx+1);
375 Double_t nEvtGen = hMCvertex->GetBinContent(ivtx+1);
376 cout<<"iBinVtx: "<<ivtx<<" nEvtsGen-> "<<nEvtGen<<" nEvtsCorr-> "<<nEvtsivtxCorr<<endl;
377 for (Int_t jeta=0; jeta<60; jeta++) {
378 if (hSPDEta_2D_bkgCorr->GetBinContent(jeta+1,ivtx+1)!=0.)
379 hnorm_bkgCorr->Fill(hBackgroundCorr->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtx);
380 if (hSPDEta_2D_bkgEffCorr->GetBinContent(jeta+1,ivtx+1)!=0.)
381 hnorm_bkgEffCorr->Fill(hTrackletRecEff->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtx);
382 if (hSPDEta_2D_bkgEffAccCorr->GetBinContent(jeta+1,ivtx+1)!=0.)
383 hnorm_bkgEffAccCorr->Fill(hBackgroundCorr->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtx);
384 if (hSPDEta_2D_fullyCorr->GetBinContent(jeta+1,ivtx+1)!=0.)
385 hnorm_fullyCorr->Fill(hBackgroundCorr->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtx);
386 if (hSPDEta_2D_triggeredEvts->GetBinContent(jeta+1,ivtx+1)!=0.)
387 hnorm_triggered->Fill(hBackgroundCorr->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtxCorr_triggered);
388 if (hSPDEta_2D_allEvts->GetBinContent(jeta+1,ivtx+1)!=0.)
389 hnorm->Fill(hBackgroundCorr->GetXaxis()->GetBinCenter(jeta+1),nEvtsivtxCorr);
390 }
391 }
392
393 for (Int_t jeta=0; jeta<60; jeta++) {
394 hnorm_bkgCorr->SetBinError(jeta+1,TMath::Sqrt(hnorm_bkgCorr->GetBinContent(jeta+1)));
395 hnorm_bkgEffCorr->SetBinError(jeta+1,TMath::Sqrt(hnorm_bkgEffCorr->GetBinContent(jeta+1)));
396 hnorm_bkgEffAccCorr->SetBinError(jeta+1,TMath::Sqrt(hnorm_bkgEffAccCorr->GetBinContent(jeta+1)));
397 hnorm_fullyCorr->SetBinError(jeta+1,TMath::Sqrt(hnorm_fullyCorr->GetBinContent(jeta+1)));
398 hnorm_triggered->SetBinError(jeta+1,TMath::Sqrt(hnorm_triggered->GetBinContent(jeta+1)));
399 hnorm->SetBinError(jeta+1,TMath::Sqrt(hnorm->GetBinContent(jeta+1)));
400 }
401
402 hdNdEta_bkgCorr->Divide(hSPDEta_2D_bkgCorr->ProjectionX("SPDEta_2D_bkgCorr_x",binVtxStart+1,binVtxStop,"e"),hnorm_bkgCorr);
403 hdNdEta_bkgEffCorr->Divide(hSPDEta_2D_bkgEffCorr->ProjectionX("SPDEta_2D_bkgEffCorr_x",binVtxStart+1,binVtxStop,"e"),hnorm_bkgEffCorr);
404 hdNdEta_bkgEffAccCorr->Divide(hSPDEta_2D_bkgEffAccCorr->ProjectionX("SPDEta_2D_bkgEffAccCorr_x",binVtxStart+1,binVtxStop,"e"),hnorm_bkgEffAccCorr);
405 hdNdEta_fullyCorr->Divide(hSPDEta_2D_fullyCorr->ProjectionX("SPDEta_2D_fullyCorr_x",binVtxStart+1,binVtxStop,"e"),hnorm_fullyCorr);
406 hdNdEta_triggeredEvts->Divide(hSPDEta_2D_triggeredEvts->ProjectionX("SPDEta_2D_triggeredEvts_x",binVtxStart+1,binVtxStop,"e"),hnorm_triggered);
407 hdNdEta_allEvts->Divide(hSPDEta_2D_allEvts->ProjectionX("SPDEta_2D_allEvts_x",binVtxStart+1,binVtxStop,"e"),hnorm);
408
409 hdNdEta_bkgCorr->Scale(10.);
410 hdNdEta_bkgEffCorr->Scale(10.);
411 hdNdEta_bkgEffAccCorr->Scale(10.);
412 hdNdEta_fullyCorr->Scale(10.);
413 hdNdEta_triggeredEvts->Scale(10.);
414 hdNdEta_allEvts->Scale(10.);
415
416 TH2F* hRatiodNdEta2D=new TH2F("ratiodNdEta2D","",60,-3,3,20,-20,20);
417// hRatiodNdEta2D->Divide(hSPDEta_2D_fullyCorr,Eta_RV,1.,1.);
418 hRatiodNdEta2D->Divide(hSPDEta_2D_allEvts,Eta,1.,1.);
419
420 TH1F* hRationorm=new TH1F("rationorm","",60,-3,3);
421 hRationorm->Add(hnorm);
422 hRationorm->Scale(1./nEvtsTot);
423
424 TH1D *hdNdEta_test = new TH1D("dNdEta_test","Pseudorapidity ",60,-3,3);
425 TH2F *hMask_allEvts = new TH2F("Mask_allEvts","",60,-3,3,20,-20,20);
426 hMask_allEvts->Sumw2();
427
428 TH2F *Eta_test = new TH2F("Eta_test","",60,-3.,3.,20,-20.,20.);
429 if (kPrimariesTest) {
430 for (Int_t kvtx=0; kvtx<20; kvtx++) {
431 for (Int_t jeta=0; jeta<60; jeta++) {
432 if (hSPDEta_2D_allEvts->GetBinContent(jeta+1,kvtx+1)!=0) {
433 hMask_allEvts->SetBinContent(jeta+1,kvtx+1,1);
434 }
435 }
436 }
437 //dN/dEta gen test
438 Eta_test->Multiply(Eta,hMask_allEvts,1.,1.);
439 hdNdEta_test->Divide(Eta_test->ProjectionX("eta_x",0,-1,"e"),hnorm,1.,1.);
440 hdNdEta_test->Scale(10.);
441 }
442
443 //Ratios generated/corrected distributions
444 TH1F* hRatiodNdEta=new TH1F("ratiodNdEta","",60,-3,3);
445 hRatiodNdEta->Sumw2();
446
447 hRatiodNdEta->SetMarkerColor(2);
448 hRatiodNdEta->SetMarkerStyle(22);
449 hRatiodNdEta->GetXaxis()->SetTitle("Pseudorapidity #eta");
450 hRatiodNdEta->GetYaxis()->SetTitle("Generated/Corrected");
451
452 TH1F* hRatiodNdEta_vertex=new TH1F("ratiodNdEta_vertex","",60,-3,3);
453 TH1F* hRatiodNdEta_triggered=new TH1F("ratiodNdEta_triggered","",60,-3,3);
454
455
456 for (Int_t i=0;i<60;i++) hdNdEta->SetBinError(i+1,0.);
457 for (Int_t i=0;i<60;i++) hdNdEta_RV->SetBinError(i+1,0.);
458
459 if (kPrimariesTest) {
460 hRatiodNdEta->Divide(hdNdEta_test,hdNdEta_allEvts,1.,1.);
461
462 } else {
463 hRatiodNdEta->Divide(hdNdEta,hdNdEta_allEvts,1.,1.);
464 hRatiodNdEta_vertex->Divide(hdNdEta_RV,hdNdEta_fullyCorr,1.,1.);
465 hRatiodNdEta_triggered->Divide(hdNdEta_Trigg,hdNdEta_triggeredEvts,1.,1.);
466 }
467
468 // Draw dN/dEta
469 hdNdEta_raw->SetLineColor(kGreen);
470 hdNdEta_raw->SetLineWidth(3);
471
472 hdNdEta_RV->SetLineWidth(3);
473 hdNdEta_RV->SetMinimum(0.);
474 hdNdEta_RV->SetMaximum(15.);
475 hdNdEta_RV->SetLineColor(kRed);
476// hdNdEta_RV->Draw("histo");
477 hdNdEta->SetLineWidth(3);
478 hdNdEta->SetMinimum(0.);
479 hdNdEta->SetMaximum(15.);
480 hdNdEta->Draw("histo");
481// hdNdEta_RV->Draw("histo,same");
482
483 hdNdEta_raw->Draw("histo,same");
484
485 hdNdEta_bkgCorr->SetMarkerStyle(21);
486 hdNdEta_bkgCorr->SetMarkerColor(kBlue);
487
488 hdNdEta_bkgCorr->Draw("same,e");
489
490 hdNdEta_bkgEffCorr->SetMarkerStyle(22);
491 hdNdEta_bkgEffCorr->SetMarkerColor(kGreen);
492
493 hdNdEta_bkgEffAccCorr->SetMarkerStyle(23);
494 hdNdEta_bkgEffAccCorr->SetMarkerColor(kRed);
495
496 hdNdEta_bkgEffCorr->Draw("p,same,e");
497 hdNdEta_bkgEffAccCorr->Draw("p,same,e");
498 hdNdEta_fullyCorr->SetMarkerStyle(20);
499 hdNdEta_fullyCorr->SetMarkerColor(kBlue);
500 hdNdEta_fullyCorr->Draw("p,same,e");
501
502 hdNdEta_triggeredEvts->SetMarkerStyle(21);
503 hdNdEta_triggeredEvts->SetMarkerColor(kOrange);
504 hdNdEta_triggeredEvts->Draw("p,same,e");
505
506 hdNdEta_allEvts->SetMarkerStyle(21);
507 hdNdEta_allEvts->SetMarkerColor(kRed);
508 hdNdEta_allEvts->Draw("p,same,e");
509
510 TLegend *leg1 = new TLegend(0.3,0.7,0.7,0.9,NULL,"brNDC");
511 leg1->SetTextFont(62);
512 leg1->SetTextSize(0.0304569);
513 TLegendEntry *entry1=leg1->AddEntry(hdNdEta_raw,"Reconstructed","l");
514// entry1=leg1->AddEntry(hdNdEta_RV,"Generated","l");
515 entry1=leg1->AddEntry(hdNdEta,"Generated","l");
516 entry1=leg1->AddEntry(hdNdEta_bkgCorr,"Bkg corrected","p");
517 entry1=leg1->AddEntry(hdNdEta_bkgEffCorr,"Bkg+AlgEff+SPDEff corrected","p");
518 entry1=leg1->AddEntry(hdNdEta_bkgEffAccCorr,"Bkg+AlgEff+SPDAccEff corrected","p");
519 entry1=leg1->AddEntry(hdNdEta_fullyCorr,"Corrected (triggered events with vertex)","p");
520 entry1=leg1->AddEntry(hdNdEta_triggeredEvts,"Corrected (triggered events)","p");
521 entry1=leg1->AddEntry(hdNdEta_allEvts,"Corrected (all events)","p");
522
523 leg1->Draw();
524
525// hRatiodNdEta->Draw();
526
527 // Write histos
528
529 TFile *fout= new TFile("SPDdNdEtaCorr.root","RECREATE");
530
531 hAccCorr->Write();
532 hAccCorrPb->Write();
533
534 hBackgroundCorr->Write();
535 hTrackletRecEff->Write();
536 hDecPartLay2Corr->Write();
537
538 hSPDEta_2Draw_accBin->Write();
539 hSPDEta_2Draw->Write();
540
541 hSPDvertex->Write();
542 hSPDvertexCorr_y->Write();
543 hSPDEta_2D_bkgCorr->Write();
544 hSPDEta_2D_bkgEffCorr->Write();
545 hSPDEta_2D_bkgEffAccCorr->Write();
546 hSPDEta_2D_fullyCorr->Write();
547 hSPDEta_2D_triggeredEvts->Write();
548 hSPDEta_2D_allEvts->Write();
549
550 hnorm_bkgCorr->Write();
551 hnorm_bkgEffCorr->Write();
552 hnorm_bkgEffAccCorr->Write();
553 hnorm_triggered->Write();
554 hnorm->Write();
555 hnorm_gen->Write();
556
557 hVertexCorr->Write();
558 hVertexTrigCorr->Write();
559 hTriggerCorr->Write();
560 hVertexCorrEta->Write();
561 hTriggerCorrEta->Write();
562 hVertexTrigCorrEta->Write();
563
564 fHistZdistribTrigVtxEvt->Write();
565 fHistEvtCorrFactor->Write();
566
567 hESDVertexvsNtracklets->Write();
568 hESDVertexvsNtracklets_trigEvts->Write();
569 hESDVertexvsNtracklets_Corr->Write();
570
571 hdNdEta->Write();
572 hdNdEta_RV->Write();
573 hdNdEta_bkgCorr->Write();
574 hdNdEta_bkgEffCorr->Write();
575 hdNdEta_bkgEffAccCorr->Write();
576 hdNdEta_fullyCorr->Write();
577 hdNdEta_triggeredEvts->Write();
578 hdNdEta_allEvts->Write();
579 hdNdEta_raw->Write();
580
581 hAccErrorsvsAccTr->Write();
582
583 hRatiodNdEta->Write();
584 hRatiodNdEta_triggered->Write();
585 hRatiodNdEta_vertex->Write();
586 hRatiodNdEta2D->Write();
587 hRationorm->Write();
588
589 Eta_test->Write();
590 Eta->Write();
591
592 hprod->Write();
593
594 fout->Write();
595 fout->Close();
596
597}