]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/QATasks/post/multistrangeQA.C
Removing obsolete task and modifying Phi PostProcess macro
[u/mrichter/AliRoot.git] / PWGLF / QATasks / post / multistrangeQA.C
CommitLineData
c683985a 1//////////////////////////////////////////////////
2//
3// This macro was written by Domenico Colella (domenico.colella@cern.ch).
4// 12 November 2013
5//
6// ------ Arguments
7// -- icasType = 0) Xi- 1) Xi+ 2) Omega- 3) Omega+
8// -- collidingsystem = 0) PbPb 1) pp
9// -- fileDir = "Input file directory"
10// -- filein = "Input file name"
11//
12// ------ QATask output content
13// The output produced by the QATask is a CFContainer with 4 steps and 21 variables.
14// The meaning of each variable within the container are listed here:
15// -- 0 = Max DCA Cascade Daughters pp: 2.0 PbPb: 0.3
16// -- 1 = Min DCA Bach To PV pp: 0.01 PbPb: 0.03
17// -- 2 = Min Cascade Cosine Of PA pp: 0.98 PbPb: 0.999
18// -- 3 = Min Cascade Radius Fid. Vol. pp: 0.2 PbPb: 0.9
19// -- 4 = Window Invariant Mass Lambda pp: 0.008 PbPb: 0.0008
20// -- 5 = Max DCA V0 Daughters pp: 1.5 PbPb: 1.0
21// -- 6 = Min V0 Cosine Of PA To PV pp: pT dep. PbPb: 0.98
22// -- 7 = Min V0 Radius Fid. Vol. pp: 0.2 PbPb: 0.9
23// -- 8 = Min DCA V0 To PV pp: 0.01 PbPb: 0.05
24// -- 9 = Min DCA Pos To PV pp: 0.05 PbPb: 0.1
25// -- 10 = Min DCA Neg To PV pp: 0.05 PbPb: 0.1
26// -- 11 = Invariant Mass distribution for Xi
27// -- 12 = Invariant Mass distribution for Omega
28// -- 13 = Transverse Momentum distribution
29// -- 14 = Rapidity distribution for Xi
30// -- 15 = Rapidity distribution for Omega
31// -- 16 = Proper length distribution for the cascade
32// -- 17 = Proper length distribution for the V0
33// -- 18 = Min V0 Cosine Of PA To Xi Vertex pp: pT dep. PbPb: pT dep.
34// -- 19 = Centrality
35// -- 20 = ESD track multiplicity
36// The last two variables are empty in the case proton-proton collisions.
37//
38// ------ Present Macro Check
39// With this macro are produced the plots of the distributions for the topological
40// variables used during the reconstruction of the cascades and defined in the two
41// classes AliCascadeVertexer.cxx and AliV0vertexer.cxx contained in /STEER/ESD/
42// folder in Aliroot.
43//
44// -- First Canvas: DCA cascade daughters, Bachelor IP to PV, Cascade cosine of PA
45// Cascade radius of fiducial volume, Invariant mass Lambda,
46// DCA V0 daughters.
47// -- Second Canvas: V0 cosine of PA to PV, Min V0 Radius Fid. Vol., Min DCA V0 To PV
48// Min DCA Pos To PV, Min DCA Neg To PV, V0 cosine of PA to XiV
49//
50// In this plots, in correspondence to the min/max cut value adopted in the reconstruction
51// a line is drawn, to check if there is same problem in the cuts definition.
52//
53// Also, other specific distribution for the selected cascades are ploted as: the
54// invariant mass distribution, the transverse momentum distribution, the rapidity
55// distribution, proper length distribution for the cascade and the v0.
56//
57// -- Third Canvas: InvMass, Transverse momentum, Cascade proper length, V0 proper length
58//
59// In the end, only for thr PbPb collision the centrality distribution and the
60// track multiplicity distribution are sored.
61//
62// -- Fourth Canvas: Centrality, track multiplicity
63//
64//
65//////////////////////////////////////////////////////
66
67
68
69
70class AliCFContainer;
71
72void multistrangeQA(Int_t icasType = 0, // 0) Xi- 1) Xi+ 2) Omega- 3) Omega+
73 Int_t collidingsystem = 0, // 0) PbPb 1) pp
74 Char_t *fileDir = "./", // Input file directory
75 Char_t *filein = "AnalysisResults_AOD.root" // Input file name
76 ) {
77
78
79 /////////////
80 gStyle->SetOptStat(1110);
81 gStyle->SetOptStat(kFALSE);
82 gStyle->SetOptTitle(kFALSE);
83 gStyle->SetFrameLineWidth(2.5);
84 gStyle->SetCanvasColor(0);
85 gStyle->SetPadColor(0);
86 gStyle->SetHistLineWidth(2.5);
87 gStyle->SetLabelSize(0.05, "x");
88 gStyle->SetLabelSize(0.05, "y");
89 gStyle->SetTitleSize(0.05, "x");
90 gStyle->SetTitleSize(0.05, "y");
91 gStyle->SetTitleOffset(1.1, "x");
92 gStyle->SetPadBottomMargin(0.14);
93 gSystem->Load("libANALYSIS.so");
94 gSystem->Load("libANALYSISalice.so");
95 gSystem->Load("libCORRFW.so");
96
97
98
99 TFile *f1 = new TFile(Form("%s/%s",fileDir,filein));
100 AliCFContainer *cf = (AliCFContainer*) (f1->Get("PWGLFStrangeness.outputCheckCascade/fCFContCascadeCuts"));
101
102
103 //DEEFINE TEXT
104 TLatex* t1 = new TLatex(0.6,0.55,"#color[3]{OK!!}");
105 t1->SetTextSize(0.1);
106 t1->SetNDC();
107 TLatex* t2 = new TLatex(0.6,0.55,"#color[2]{NOT OK!!}");
108 t2->SetTextSize(0.1);
109 t2->SetNDC();
110 t2->SetTextColor(2);
111 TLatex* tcasc;
112 if (icasType == 0) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{-}}");
113 else if (icasType == 1) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{+}}");
114 else if (icasType == 2) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{-}}");
115 else if (icasType == 3) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{+}}");
116 tcasc->SetTextSize(0.2);
117 tcasc->SetNDC();
118 tcasc->SetTextColor(2);
119 TLatex* tpdgmass;
120 if (icasType == 0) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
121 else if (icasType == 1) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
122 else if (icasType == 2) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
123 else if (icasType == 3) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
124 tpdgmass->SetTextSize(0.07);
125 tpdgmass->SetNDC();
126 tpdgmass->SetTextColor(2);
127
128 //DEFINE 1st CANVAS AND DRAW PLOTS
129 TCanvas *c1 = new TCanvas("c1","",1200,800);
130 c1->Divide(2,3);
131 //Pad 1: DCA cascade daughters
132 c1->cd(1);
133 gPad->SetLogy();
134 TH1D *hvar0 = cf->ShowProjection(0,icasType);
135 hvar0->Draw("histo");
136 Double_t x0;
137 if (collidingsystem == 0) x0 = 0.3;
138 else if (collidingsystem == 1) x0 = 2.0;
139 TLine *line0 = new TLine(x0,0.,x0,hvar0->GetBinContent(hvar0->GetMaximumBin()));
140 line0->SetLineColor(kRed);
141 line0->SetLineStyle(9);
142 line0->SetLineWidth(2.0);
143 line0->Draw("same");
144 Bool_t check_0 = checkOverTheLimit(hvar0, x0);
145 if (check_0) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
146 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
147 tcasc->Draw();
148 //Pad 2: Bachelor IP to PV
149 c1->cd(2);
150 gPad->SetLogy();
151 TH1D *hvar1 = cf->ShowProjection(1,icasType);
152 hvar1->GetXaxis()->SetRangeUser(0.,0.24);
153 hvar1->Draw("histo");
154 Double_t x1;
155 if (collidingsystem == 0) x1 = 0.03;
156 else if (collidingsystem == 1) x1 = 0.01;
157 TLine *line1 = new TLine(x1,0.,x1,hvar1->GetBinContent(hvar1->GetMaximumBin()));
158 line1->SetLineColor(kRed);
159 line1->SetLineStyle(9);
160 line1->SetLineWidth(2.0);
161 line1->Draw("same");
162 Bool_t check_1 = checkUnderTheLimit(hvar1, x1);
163 if (check_1) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
164 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
165 //Pad 3: Cascade cosine of Pointing Angle
166 c1->cd(3);
167 gPad->SetLogy();
168 TH1D *hvar2 = cf->ShowProjection(2,icasType);
169 Double_t max2 = hvar2->GetBinContent(hvar2->GetMaximumBin());
170 hvar2->GetYaxis()->SetRangeUser(0.01,max2*1.5);
171 hvar2->Draw("histo");
172 Double_t x2;
173 if (collidingsystem == 0) x2 = 0.999;
174 else if (collidingsystem == 1) x2 = 0.98;
175 TLine *line2 = new TLine(x2,0.,x2,hvar2->GetBinContent(hvar2->GetMaximumBin()));
176 line2->SetLineColor(kRed);
177 line2->SetLineStyle(9);
178 line2->SetLineWidth(2.0);
179 line2->Draw("same");
180 line1->Draw("same");
181 Bool_t check_2 = checkUnderTheLimit(hvar2, x2);
182 if (check_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
183 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
184 //Pad 4: Cascade radius of fiducial volume
185 c1->cd(4);
186 gPad->SetLogy();
187 TH1D *hvar3 = cf->ShowProjection(3,icasType);
188 hvar3->GetXaxis()->SetRangeUser(0.,3.8);
189 hvar3->Draw("histo");
190 Double_t x3;
191 if (collidingsystem == 0) x3 = 0.9;
192 else if (collidingsystem == 1) x3 = 0.2;
193 TLine *line3 = new TLine(x3,0.,x3,hvar3->GetBinContent(hvar3->GetMaximumBin()));
194 line3->SetLineColor(kRed);
195 line3->SetLineStyle(9);
196 line3->SetLineWidth(2.0);
197 line3->Draw("same");
198 Bool_t check_3 = checkUnderTheLimit(hvar3, x3);
199 if (check_3) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
200 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
201 //Pad 5: Invariant mass Lambda
202 c1->cd(5);
203 TH1D *hvar4 = cf->ShowProjection(4,icasType);
204 hvar4->Draw("histo");
205 Double_t x41 = 1.116 + 0.008;
206 TLine *line41 = new TLine(x41,0.,x41,hvar4->GetBinContent(hvar4->GetMaximumBin()));
207 line41->SetLineColor(kRed);
208 line41->SetLineStyle(9);
209 line41->SetLineWidth(2.0);
210 line41->Draw("same");
211 Double_t x42 = 1.115 - 0.008;
212 TLine *line42 = new TLine(x42,0.,x42,hvar4->GetBinContent(hvar4->GetMaximumBin()));
213 line42->SetLineColor(kRed);
214 line42->SetLineStyle(9);
215 line42->SetLineWidth(2.0);
216 line42->Draw("same");
217 Bool_t check_4_1 = checkUnderTheLimit(hvar3, x3);
218 Bool_t check_4_2 = checkOverTheLimit(hvar0, x0);
219 if (check_4_1 && check_4_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
220 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
221 //Pad 6: DCA V0 daughters
222 c1->cd(6);
223 gPad->SetLogy();
224 TH1D *hvar5 = cf->ShowProjection(5,icasType);
225 hvar5->Draw("histo");
226 Double_t x5;
227 if (collidingsystem == 0) x5 = 1.0;
228 else if (collidingsystem == 1) x5 = 1.5;
229 TLine *line5 = new TLine(x5,0.,x5,hvar5->GetBinContent(hvar5->GetMaximumBin()));
230 line5->SetLineColor(kRed);
231 line5->SetLineStyle(9);
232 line5->SetLineWidth(2.0);
233 line5->Draw("same");
234 Bool_t check_5 = checkOverTheLimit(hvar5, x5);
235 if (check_5) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
236 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
237 c1->SaveAs("MultistrangeQA.pdf(");
238
239
240 //DEFINE 2st CANVAS AND DRAW PLOTS
241 TCanvas *c2 = new TCanvas("c2","",1200,800);
242 c2->Divide(2,3);
243 //Pad 1: V0 cosine of Pointing Angle to PV
244 c2->cd(1);
245 gPad->SetLogy();
246 TH1D *hvar6 = cf->ShowProjection(6,icasType);
247 Double_t max6 = hvar6->GetBinContent(hvar6->GetMaximumBin());
248 hvar6->GetYaxis()->SetRangeUser(0.01,max6*1.5);
249 hvar6->Draw("histo");
250 //Pad 2: Min V0 Radius Fid. Vol.
251 c2->cd(2);
252 gPad->SetLogy();
253 TH1D *hvar7 = cf->ShowProjection(7,icasType);
254 hvar7->GetXaxis()->SetRangeUser(0.,3.0);
255 hvar7->Draw("histo");
256 Double_t x7;
257 if (collidingsystem == 0) x7 = 0.9;
258 else if (collidingsystem == 1) x7 = 0.2;
259 TLine *line7 = new TLine(x7,0.,x7,hvar7->GetBinContent(hvar7->GetMaximumBin()));
260 line7->SetLineColor(kRed);
261 line7->SetLineStyle(9);
262 line7->SetLineWidth(2.0);
263 line7->Draw("same");
264 Bool_t check_7 = checkUnderTheLimit(hvar7, x7);
265 if (check_7) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
266 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
267 //Pad3: Min DCA V0 To PV
268 c2->cd(3);
269 gPad->SetLogy();
270 TH1D *hvar8 = cf->ShowProjection(8,icasType);
271 hvar8->GetXaxis()->SetRangeUser(0.,0.3);
272 hvar8->Draw("histo");
273 Double_t x8;
274 if (collidingsystem == 0) x8 = 0.05;
275 else if (collidingsystem == 1) x8 = 0.01;
276 TLine *line8 = new TLine(x8,0.,x8,hvar8->GetBinContent(hvar8->GetMaximumBin()));
277 line8->SetLineColor(kRed);
278 line8->SetLineStyle(9);
279 line8->SetLineWidth(2.0);
280 line8->Draw("same");
281 Bool_t check_8 = checkUnderTheLimit(hvar8, x8);
282 if (check_8) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
283 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
284 //Pad 4: Min DCA Pos To PV
285 c2->cd(4);
286 gPad->SetLogy();
287 TH1D *hvar9 = cf->ShowProjection(9,icasType);
288 hvar9->GetXaxis()->SetRangeUser(0.,0.2);
289 hvar9->Draw("histo");
290 Double_t x9;
291 if (collidingsystem == 0) x9 = 0.1;
292 else if (collidingsystem == 1) x9 = 0.05;
293 TLine *line9 = new TLine(x9,0.,x9,hvar9->GetBinContent(hvar9->GetMaximumBin()));
294 line9->SetLineColor(kRed);
295 line9->SetLineStyle(9);
296 line9->SetLineWidth(2.0);
297 line9->Draw("same");
298 Bool_t check_9 = checkUnderTheLimit(hvar9, x9);
299 if (check_9) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
300 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
301 //Pad 5: Min DCA Neg To PV
302 c2->cd(5);
303 gPad->SetLogy();
304 TH1D *hvar10 = cf->ShowProjection(10,icasType);
305 hvar10->GetXaxis()->SetRangeUser(0.,0.2);
306 hvar10->Draw("histo");
307 Double_t x10;
308 if (collidingsystem == 0) x10 = 0.1;
309 else if (collidingsystem == 1) x10 = 0.05;
310 TLine *line10 = new TLine(x10,0.,x10,hvar10->GetBinContent(hvar10->GetMaximumBin()));
311 line10->SetLineColor(kRed);
312 line10->SetLineStyle(9);
313 line10->SetLineWidth(2.0);
314 line10->Draw("same");
315 Bool_t check_10 = checkUnderTheLimit(hvar10, x10);
316 if (check_10) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
317 else { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
318 //Pad 6: V0 cosine of Pointing Angle to XiV
319 c2->cd(6);
320 gPad->SetLogy();
321 TH1D *hvar20 = cf->ShowProjection(18,icasType);
322 Double_t max20 = hvar20->GetBinContent(hvar20->GetMaximumBin());
323 hvar20->GetYaxis()->SetRangeUser(0.01,max20*1.5);
324 hvar20->Draw("histo");
325 c2->SaveAs("MultistrangeQA.pdf");
326
327 //DEFINE 3st CANVAS AND DRAW PLOTS
328 TCanvas *c3 = new TCanvas("c3","",1200,800);
329 c3->Divide(2,3);
330 //Pad 1: InvMass
331 c3->cd(1);
332 TH1D *hvar12 = cf->ShowProjection(11+icasType/2,icasType);
333 hvar12->Draw("histo");
334 tpdgmass->Draw();
335 TLine *linemass;
336 if (icasType == 0) linemass = new TLine(1.32171,0.,1.32171,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
337 else if (icasType == 1) linemass = new TLine(1.32171,0.,1.32171,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
338 else if (icasType == 2) linemass = new TLine(1.67245,0.,1.67245,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
339 else if (icasType == 3) linemass = new TLine(1.67245,0.,1.67245,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
340 linemass->SetLineColor(kRed);
341 linemass->SetLineStyle(1);
342 linemass->SetLineWidth(2.0);
343 linemass->Draw("same");
344 //Pad 2: Transverse momentum
345 c3->cd(2);
346 TH1D *hvar13 = cf->ShowProjection(13,icasType);
347 hvar13->Draw("histo");
348 //Pad 3: Y
349 c3->cd(3);
350 TH1D *hvar14 = cf->ShowProjection(14+icasType/2,icasType);
351 hvar14->Draw("histo");
352 //Pad 4: Cascade proper length
353 c3->cd(4);
354 TH1D *hvar18;
355 hvar18 = cf->ShowProjection(16,icasType);
356 hvar18->GetXaxis()->SetRangeUser(0.,90.);
357 hvar18->Draw("histo");
358 //Pad 5: V0 proper length
359 c3->cd(5);
360 TH1D *hvar19;
361 hvar19 = cf->ShowProjection(17,icasType);
362 hvar19->GetXaxis()->SetRangeUser(0.,90.);
363 hvar19->Draw("histo");
364 //Pad 6
365 // empty
366 if (collidingsystem == 1) c3->SaveAs("MultistrangeQA.pdf)");
367 else if (collidingsystem == 0) c3->SaveAs("MultistrangeQA.pdf");
368
369
370 //DEFINE 4st CANVAS AND DRAW PLOTS
371 TCanvas *c4 = new TCanvas("c4","",600,400);
372 c4->Divide(2,1);
373 //Pad1: invariant mass fit
374 c4->cd(1);
375 TH1D *hvar18 = cf->ShowProjection(11+icasType/2,icasType);
376 hvar18->Draw("histo");
377 // - SOME PARAMETER VALUE
378 Bool_t kfitgauss = kFALSE;
379 Bool_t kfitleft = kFALSE;
380 Bool_t kfitright = kFALSE;
381 Int_t ptbinNarrowY = 0;
382 if (icasType < 2) ptbinNarrowY = 10; // 6;
383 else ptbinNarrowY = 3; // 2;
384 // - SOME DEFINITIONS
385 Float_t lowlimmass;
386 Float_t uplimmass;
387 Float_t lowgausslim;
388 Float_t upgausslim;
389 if (icasType==0||icasType==1) {
390 lowlimmass=1.30;
391 uplimmass=1.34;
392 lowgausslim=1.312;
393 upgausslim=1.332;
394 } else {
395 lowlimmass=1.645;
396 uplimmass=1.70;
397 lowgausslim=1.668;
398 upgausslim=1.678;
399 }
400 TF1* fitinvmass = new TF1("fitinvmass","gaus(0)+pol2(3)",lowlimmass,uplimmass);
401 fitinvmass->SetParName(0, "cnstntG");
402 fitinvmass->SetParName(1, "meanG");
403 fitinvmass->SetParName(2, "sigmaG");
404 fitinvmass->SetParLimits(0,0.,500000.);
405 if (icasType==0||icasType==1) {
406 fitinvmass->SetParameter(1, 1.32171);
407 fitinvmass->SetParLimits(1, 1.31,1.33);
408 fitinvmass->SetParLimits(2,0.001,0.005);
409 } else {
410 fitinvmass->SetParameter(1, 1.67245);
411 fitinvmass->SetParLimits(1, 1.664,1.68);
412 fitinvmass->SetParLimits(2,0.0008,0.006);
413 }
414 hvar18->Fit("fitinvmass","rimeN");
415 fitinvmass->SetLineColor(kRed);
416 fitinvmass->Draw("same");
417 Float_t meanGauss = fitinvmass->GetParameter(1);
418 Float_t sigmaGauss = fitinvmass->GetParameter(2);
419 cout<<"Mean: "<<meanGauss<<endl;
420 cout<<"Sigma: "<<sigmaGauss<<endl;
421 //Pad2: Text
422 c4->cd(2);
423 Float_t refwidth = 0.002;
424 TPaveText *pave1 = new TPaveText(0.05,0.3,0.95,0.5);
425 pave1->SetFillColor(0);
426 pave1->SetTextSize(0.04);
427 pave1->SetTextAlign(12);
428 if (icasType < 2) pave1->AddText("PDG mass: 1.32171 GeV/c^{2}");
429 else pave1->AddText("PDG mass: 1.67245 GeV/c^{2}");
430 pave1->AddText(Form("#color[1]{Mass form Fit: %.5f #pm %.5f GeV/c^{2}}",meanGauss,sigmaGauss));
431 if (sigmaGauss > refwidth - 0.0003 && sigmaGauss < refwidth + 0.0003) pave1->AddText("#color[3]{OK!! The width is compatible with standard.}");
432 else pave1->AddText("#color[2]{NOT OK!! Problem.}");
433 pave1->Draw();
434 cout<<" "<<refwidth - 0.0003<<"<"<<sigmaGauss<<"<"<<refwidth + 0.0003<<endl;
435
436 //DEFINE 5st CANVAS AND DRAW PLOTS
437 if (collidingsystem == 0) {
438 TCanvas *c5 = new TCanvas("c5","",1200,270);
439 c5->Divide(2,1);
440 //Pad 1: centrality
441 c5->cd(1);
442 TH1D *hvar16 = cf->ShowProjection(19,icasType);
443 hvar16->Draw("histo");
444 //Pad 2: track multiplicity
445 c5->cd(2);
446 TH1D *hvar17 = cf->ShowProjection(20,icasType);
447 hvar17->Draw("histo");
448 c5->SaveAs("MultistrangeQA.pdf)");
449 }
450
451
452}
453
454
455
456
457Bool_t checkUnderTheLimit(TH1D *lHist, Double_t limit) {
458
459 Int_t binlimit = lHist->FindBin(limit);
460 Bool_t checkOk = kTRUE;
461 for (Int_t i = 1; i < binlimit; i++) {
462 Int_t content = 0;
463 content = lHist->GetBinContent(i);
464 if (content != 0) checkOk = kFALSE;
465 //cout<<"Content bin "<<i<<": "<<content<<endl;
466 }
467 return checkOk;
468
469}
470
471
472Bool_t checkOverTheLimit(TH1D *lHist, Double_t limit) {
473
474 Int_t binlimit = lHist->FindBin(limit);
475 Int_t lastbin = lHist->GetNbinsX();
476 Bool_t checkOk = kTRUE;
477 for (Int_t i = binlimit; i < lastbin+1; i++) {
478 Int_t content = 0;
479 content = lHist->GetBinContent(i);
480 if (content != 0) checkOk = kFALSE;
481 //cout<<"Content bin "<<i<<": "<<content<<endl;
482 }
483 return checkOk;
484
485}
486
487