]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/analysisQA/processMultistrangeQA.C
.so cleanup: removed from gSystem->Load()
[u/mrichter/AliRoot.git] / PWGPP / analysisQA / processMultistrangeQA.C
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
70 class AliCFContainer;
71
72 void processMultistrangeQA(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.root",   // Input file name
76                            TString suffix          = "eps",
77                            const char *outfile     = "output.root"             // output
78                            ) {
79
80
81       /////////////
82       gStyle->SetOptStat(1110);
83       gStyle->SetOptStat(kFALSE);
84       gStyle->SetOptTitle(kFALSE);
85       gStyle->SetFrameLineWidth(2.5);
86       gStyle->SetCanvasColor(0);
87       gStyle->SetPadColor(0);
88       gStyle->SetHistLineWidth(2.5);
89       gStyle->SetLabelSize(0.05, "x");
90       gStyle->SetLabelSize(0.05, "y");
91       gStyle->SetTitleSize(0.05, "x");
92       gStyle->SetTitleSize(0.05, "y");
93       gStyle->SetTitleOffset(1.1, "x");
94       gStyle->SetPadBottomMargin(0.14);
95       gSystem->Load("libANALYSIS");
96       gSystem->Load("libANALYSISalice");
97       gSystem->Load("libCORRFW");
98
99  
100
101      TFile *f1 = new TFile(Form("%s/%s",fileDir,filein));
102      AliCFContainer *cf = (AliCFContainer*) (f1->Get("PWGLFStrangeness.outputCheckCascade/fCFContCascadeCuts"));
103      if(!cf) {
104        Printf("FATAL: PWGLFStrangeness.outputCheckCascade/fCFContCascadeCuts ===> Not Available");
105        Printf("Exiting processMultistrangeQA");
106        return;
107      }
108
109
110      //DEEFINE TEXT
111      TLatex* t1 = new TLatex(0.6,0.55,"#color[3]{OK!!}");
112      t1->SetTextSize(0.1);
113      t1->SetNDC();
114      TLatex* t2 = new TLatex(0.6,0.55,"#color[2]{NOT OK!!}");
115      t2->SetTextSize(0.1);
116      t2->SetNDC();
117      t2->SetTextColor(2);
118      TLatex* tcasc;
119      if      (icasType == 0) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{-}}");
120      else if (icasType == 1) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{+}}");
121      else if (icasType == 2) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{-}}");
122      else if (icasType == 3) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{+}}");
123      tcasc->SetTextSize(0.2);
124      tcasc->SetNDC();
125      tcasc->SetTextColor(2);
126      TLatex* tpdgmass;
127      if      (icasType == 0) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
128      else if (icasType == 1) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
129      else if (icasType == 2) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
130      else if (icasType == 3) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
131      tpdgmass->SetTextSize(0.07);
132      tpdgmass->SetNDC();
133      tpdgmass->SetTextColor(2);
134  
135      // Added by sjena
136      TFile *fout = TFile::Open(outfile,"UPDATE");
137      fout->ls();
138      
139      TDirectoryFile *cdd = NULL;
140      cdd = (TDirectoryFile*)fout->Get("LF");
141      if(!cdd) {
142        Printf("Warning: LF <dir> doesn't exist, creating a new one");
143        cdd = (TDirectoryFile*)fout->mkdir("LF");
144      }
145      cdd->cd();
146      cdd->ls();
147      
148
149
150      // TFile *f = TFile::Open("OutPut.root","UPDATE");
151
152
153      //DEFINE 1st CANVAS AND DRAW PLOTS
154      TCanvas *c1 = new TCanvas("c1","",1200,800);
155      c1->Divide(2,3); 
156        //Pad 1: DCA cascade daughters
157        c1->cd(1);
158        gPad->SetLogy();
159        TH1D *hvar0 = cf->ShowProjection(0,icasType);
160        hvar0->Draw("histo");
161
162        //   hvar0->SetName(Form("fig_lf_multistrange_%s", hvar0->GetName()));
163        hvar0->SetName(Form("fig_lf_multistrange_0"));
164        hvar0->Write();
165
166
167
168        Double_t x0;
169        if      (collidingsystem == 0) x0 = 0.3;
170        else if (collidingsystem == 1) x0 = 2.0;
171        TLine *line0 = new TLine(x0,0.,x0,hvar0->GetBinContent(hvar0->GetMaximumBin()));
172        line0->SetLineColor(kRed);
173        line0->SetLineStyle(9);
174        line0->SetLineWidth(2.0);
175        line0->Draw("same");
176           Bool_t check_0 = checkOverTheLimit(hvar0, x0);
177           if (check_0) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
178           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
179        tcasc->Draw();
180        //Pad 2: Bachelor IP to PV
181        c1->cd(2);
182        gPad->SetLogy();
183        TH1D *hvar1 = cf->ShowProjection(1,icasType);
184        hvar1->GetXaxis()->SetRangeUser(0.,0.24);
185        hvar1->Draw("histo");
186
187        //   hvar1->SetName(Form("fig_lf_multistrange_%s", hvar1->GetName()));
188        hvar1->SetName(Form("fig_lf_multistrange_1"));
189        hvar1->Write();
190
191
192        Double_t x1;
193        if      (collidingsystem == 0) x1 = 0.03;
194        else if (collidingsystem == 1) x1 = 0.01;
195        TLine *line1 = new TLine(x1,0.,x1,hvar1->GetBinContent(hvar1->GetMaximumBin()));
196        line1->SetLineColor(kRed);
197        line1->SetLineStyle(9);
198        line1->SetLineWidth(2.0);
199        line1->Draw("same");
200           Bool_t check_1 = checkUnderTheLimit(hvar1, x1);
201           if (check_1) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
202           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
203        //Pad 3: Cascade cosine of Pointing Angle
204        c1->cd(3);
205        gPad->SetLogy();
206        TH1D *hvar2 = cf->ShowProjection(2,icasType);
207        Double_t max2 = hvar2->GetBinContent(hvar2->GetMaximumBin());
208        hvar2->GetYaxis()->SetRangeUser(0.01,max2*1.5);
209        hvar2->Draw("histo");
210
211        //     hvar2->SetName(Form("fig_lf_multistrange_%s", hvar2->GetName()));
212        hvar2->SetName(Form("fig_lf_multistrange_2"));
213        hvar2->Write();
214        
215
216        Double_t x2;
217        if      (collidingsystem == 0) x2 = 0.999;
218        else if (collidingsystem == 1) x2 = 0.98;
219        TLine *line2 = new TLine(x2,0.,x2,hvar2->GetBinContent(hvar2->GetMaximumBin()));
220        line2->SetLineColor(kRed);
221        line2->SetLineStyle(9);
222        line2->SetLineWidth(2.0);
223        line2->Draw("same");
224        line1->Draw("same");
225           Bool_t check_2 = checkUnderTheLimit(hvar2, x2);
226           if (check_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
227           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
228        //Pad 4: Cascade radius of fiducial volume
229        c1->cd(4);
230        gPad->SetLogy();
231        TH1D *hvar3 = cf->ShowProjection(3,icasType);
232        hvar3->GetXaxis()->SetRangeUser(0.,3.8);
233        hvar3->Draw("histo");
234
235        //   hvar3->SetName(Form("fig_lf_multistrange_%s", hvar3->GetName()));
236        hvar3->SetName(Form("fig_lf_multistrange_3"));
237        hvar3->Write();
238
239        Double_t x3;
240        if      (collidingsystem == 0) x3 = 0.9;
241        else if (collidingsystem == 1) x3 = 0.2;
242        TLine *line3 = new TLine(x3,0.,x3,hvar3->GetBinContent(hvar3->GetMaximumBin()));
243        line3->SetLineColor(kRed);
244        line3->SetLineStyle(9);
245        line3->SetLineWidth(2.0);
246        line3->Draw("same");
247           Bool_t check_3 = checkUnderTheLimit(hvar3, x3);
248           if (check_3) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
249           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
250        //Pad 5: Invariant mass Lambda
251        c1->cd(5);
252        TH1D *hvar4 = cf->ShowProjection(4,icasType);
253        hvar4->Draw("histo");
254
255        hvar4->SetName(Form("fig_lf_multistrange_4"));
256        hvar4->Write();
257
258        Double_t x41 = 1.116 + 0.008;
259        TLine *line41 = new TLine(x41,0.,x41,hvar4->GetBinContent(hvar4->GetMaximumBin()));
260        line41->SetLineColor(kRed);
261        line41->SetLineStyle(9);
262        line41->SetLineWidth(2.0);
263        line41->Draw("same");
264        Double_t x42 = 1.115 - 0.008;
265        TLine *line42 = new TLine(x42,0.,x42,hvar4->GetBinContent(hvar4->GetMaximumBin()));
266        line42->SetLineColor(kRed);
267        line42->SetLineStyle(9);
268        line42->SetLineWidth(2.0);
269        line42->Draw("same");
270           Bool_t check_4_1 = checkUnderTheLimit(hvar3, x3);
271           Bool_t check_4_2 = checkOverTheLimit(hvar0, x0);
272           if (check_4_1 && check_4_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
273           else                        { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
274        //Pad 6: DCA V0 daughters
275        c1->cd(6);
276        gPad->SetLogy();
277        TH1D *hvar5 = cf->ShowProjection(5,icasType);
278        hvar5->Draw("histo");
279
280        //   hvar5->SetName(Form("fig_lf_multistrange_%s", hvar5->GetName()));
281        hvar5->SetName(Form("fig_lf_multistrange_5"));
282        hvar5->Write();
283
284        Double_t x5;
285        if      (collidingsystem == 0) x5 = 1.0;
286        else if (collidingsystem == 1) x5 = 1.5;
287        TLine *line5 = new TLine(x5,0.,x5,hvar5->GetBinContent(hvar5->GetMaximumBin()));
288        line5->SetLineColor(kRed);
289        line5->SetLineStyle(9);
290        line5->SetLineWidth(2.0);
291        line5->Draw("same");
292           Bool_t check_5 = checkOverTheLimit(hvar5, x5);
293           if (check_5) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
294           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
295
296           c1->SaveAs(Form("fig_lf_MultistrangeQA_1.%s",suffix.Data()));
297           
298     
299
300      //DEFINE 2st CANVAS AND DRAW PLOTS
301      TCanvas *c2 = new TCanvas("c2","",1200,800);
302      c2->Divide(2,3);
303        //Pad 1: V0 cosine of Pointing Angle to PV
304        c2->cd(1);
305        gPad->SetLogy();
306        TH1D *hvar6 = cf->ShowProjection(6,icasType);
307        Double_t max6 = hvar6->GetBinContent(hvar6->GetMaximumBin());
308        hvar6->GetYaxis()->SetRangeUser(0.01,max6*1.5);
309        hvar6->Draw("histo");
310
311        //  hvar6->SetName(Form("fig_lf_multistrange_%s", hvar6->GetName()));
312        hvar6->SetName(Form("fig_lf_multistrange_6"));
313        hvar6->Write();
314
315        //Pad 2: Min V0 Radius Fid. Vol.  
316        c2->cd(2);
317        gPad->SetLogy();
318        TH1D *hvar7 = cf->ShowProjection(7,icasType);
319        hvar7->GetXaxis()->SetRangeUser(0.,3.0);
320        hvar7->Draw("histo");
321
322        hvar7->SetName(Form("fig_lf_multistrange_7"));
323        hvar7->Write();
324
325        Double_t x7;
326        if      (collidingsystem == 0) x7 = 0.9;
327        else if (collidingsystem == 1) x7 = 0.2;
328        TLine *line7 = new TLine(x7,0.,x7,hvar7->GetBinContent(hvar7->GetMaximumBin()));
329        line7->SetLineColor(kRed);
330        line7->SetLineStyle(9);
331        line7->SetLineWidth(2.0);
332        line7->Draw("same");
333           Bool_t check_7 = checkUnderTheLimit(hvar7, x7);
334           if (check_7) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
335           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
336        //Pad3: Min DCA V0 To PV
337        c2->cd(3);
338        gPad->SetLogy();
339        TH1D *hvar8 = cf->ShowProjection(8,icasType);
340        hvar8->GetXaxis()->SetRangeUser(0.,0.3);
341        hvar8->Draw("histo");
342
343        //      hvar8->SetName(Form("fig_lf_multistrange_%s", hvar8->GetName()));
344        hvar8->SetName(Form("fig_lf_multistrange_8"));
345        hvar8->Write();
346
347        Double_t x8;
348        if      (collidingsystem == 0) x8 = 0.05;
349        else if (collidingsystem == 1) x8 = 0.01;
350        TLine *line8 = new TLine(x8,0.,x8,hvar8->GetBinContent(hvar8->GetMaximumBin()));
351        line8->SetLineColor(kRed);
352        line8->SetLineStyle(9);
353        line8->SetLineWidth(2.0);
354        line8->Draw("same");
355           Bool_t check_8 = checkUnderTheLimit(hvar8, x8);
356           if (check_8) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
357           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
358        //Pad 4: Min DCA Pos To PV
359        c2->cd(4);
360        gPad->SetLogy();
361        TH1D *hvar9 = cf->ShowProjection(9,icasType);
362        hvar9->GetXaxis()->SetRangeUser(0.,0.2);
363        hvar9->Draw("histo");
364
365        //   hvar9->SetName(Form("fig_lf_multistrange_%s", hvar9->GetName()));
366        hvar9->SetName(Form("fig_lf_multistrange_9"));
367        hvar9->Write();
368
369        Double_t x9;
370        if      (collidingsystem == 0) x9 = 0.1;
371        else if (collidingsystem == 1) x9 = 0.05;
372        TLine *line9 = new TLine(x9,0.,x9,hvar9->GetBinContent(hvar9->GetMaximumBin()));
373        line9->SetLineColor(kRed);
374        line9->SetLineStyle(9);
375        line9->SetLineWidth(2.0);
376        line9->Draw("same");
377
378           Bool_t check_9 = checkUnderTheLimit(hvar9, x9);
379           if (check_9) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
380           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
381        //Pad 5: Min DCA Neg To PV
382        c2->cd(5);
383        gPad->SetLogy();
384        TH1D *hvar10 = cf->ShowProjection(10,icasType);
385        hvar10->GetXaxis()->SetRangeUser(0.,0.2);
386        hvar10->Draw("histo");
387
388        //    hvar10->SetName(Form("fig_lf_multistrange_%s", hvar10->GetName()));
389        hvar10->SetName(Form("fig_lf_multistrange_10"));
390        hvar10->Write();
391
392
393        Double_t x10;
394        if      (collidingsystem == 0) x10 = 0.1;
395        else if (collidingsystem == 1) x10 = 0.05;
396        TLine *line10 = new TLine(x10,0.,x10,hvar10->GetBinContent(hvar10->GetMaximumBin()));
397        line10->SetLineColor(kRed);
398        line10->SetLineStyle(9);
399        line10->SetLineWidth(2.0);
400        line10->Draw("same");
401           Bool_t check_10 = checkUnderTheLimit(hvar10, x10);
402           if (check_10) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
403           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
404        //Pad 6: V0 cosine of Pointing Angle to XiV
405        c2->cd(6);
406        gPad->SetLogy();
407        TH1D *hvar20 = cf->ShowProjection(18,icasType);
408        Double_t max20 = hvar20->GetBinContent(hvar20->GetMaximumBin());
409        hvar20->GetYaxis()->SetRangeUser(0.01,max20*1.5);
410        hvar20->Draw("histo");
411
412        //  hvar20->SetName(Form("fig_lf_multistrange_%s", hvar20->GetName()));
413        hvar20->SetName(Form("fig_lf_multistrange_20"));
414        hvar20->Write();
415
416
417        c2->SaveAs(Form("fig_lf_MultistrangeQA_2.%s",suffix.Data()));
418        
419
420      //DEFINE 3st CANVAS AND DRAW PLOTS
421      TCanvas *c3 = new TCanvas("c3","",1200,800);
422      c3->Divide(2,3);
423        //Pad 1: InvMass
424        c3->cd(1);
425        TH1D *hvar12 = cf->ShowProjection(11+icasType/2,icasType);
426        hvar12->Draw("histo");
427
428        hvar12->SetName(Form("fig_lf_multistrange_12", hvar12->GetName()));
429        hvar12->Write();
430
431        tpdgmass->Draw(); 
432        TLine *linemass;
433        if      (icasType == 0) linemass = new TLine(1.32171,0.,1.32171,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
434        else if (icasType == 1) linemass = new TLine(1.32171,0.,1.32171,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
435        else if (icasType == 2) linemass = new TLine(1.67245,0.,1.67245,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
436        else if (icasType == 3) linemass = new TLine(1.67245,0.,1.67245,0.5*hvar12->GetBinContent(hvar12->GetMaximumBin()));
437        linemass->SetLineColor(kRed);
438        linemass->SetLineStyle(1);
439        linemass->SetLineWidth(2.0);
440        linemass->Draw("same");
441        //Pad 2: Transverse momentum
442        c3->cd(2);
443        TH1D *hvar13 = cf->ShowProjection(13,icasType);
444        hvar13->Draw("histo");
445
446        //  hvar13->SetName(Form("fig_lf_multistrange_%s", hvar13->GetName()));
447        hvar13->SetName(Form("fig_lf_multistrange_13"));
448        hvar13->Write();
449        //Pad 3: Y
450        c3->cd(3);
451        TH1D *hvar14 = cf->ShowProjection(14+icasType/2,icasType);
452        hvar14->Draw("histo");
453
454        //  hvar14->SetName(Form("fig_lf_multistrange_%s", hvar14->GetName()));
455        hvar14->SetName(Form("fig_lf_multistrange_14"));
456        hvar14->Write();
457
458        //Pad 4: Cascade proper length
459        c3->cd(4);
460        TH1D *hvar18;
461        hvar18 = cf->ShowProjection(16,icasType);
462        hvar18->GetXaxis()->SetRangeUser(0.,90.);
463        hvar18->Draw("histo");
464
465        //   hvar18->SetName(Form("fig_lf_multistrange_%s", hvar18->GetName()));
466        hvar18->SetName(Form("fig_lf_multistrange_18", hvar18->GetName()));
467        hvar18->Write();
468
469        //Pad 5: V0 proper length 
470        c3->cd(5);
471        TH1D *hvar19;
472        hvar19 = cf->ShowProjection(17,icasType);
473        hvar19->GetXaxis()->SetRangeUser(0.,90.);
474        hvar19->Draw("histo");
475        
476        hvar19->SetName(Form("fig_lf_multistrange_19", hvar19->GetName()));
477        hvar19->Write();
478        
479
480       //Pad 6
481       // empty 
482        if      (collidingsystem == 1) { c3->SaveAs(Form("fig_lf_MultistrangeQA_3.%s",suffix.Data())); }
483        else if (collidingsystem == 0) { c3->SaveAs(Form("fig_lf_MultistrangeQA_3.%s",suffix.Data())); }
484
485     
486      //DEFINE 4st CANVAS AND DRAW PLOTS
487     TCanvas *c4 = new TCanvas("c4","",600,400);
488     c4->Divide(2,1);
489       //Pad1: invariant mass fit
490       c4->cd(1);
491       TH1D *hvar18 = cf->ShowProjection(11+icasType/2,icasType);
492       hvar18->Draw("histo");
493
494       hvar18->SetName(Form("fig_lf_multistrange_18_1", hvar18->GetName()));
495       hvar18->Write();
496
497        // - SOME PARAMETER VALUE
498        Bool_t kfitgauss = kFALSE;
499        Bool_t kfitleft  = kFALSE;
500        Bool_t kfitright = kFALSE;
501        Int_t  ptbinNarrowY = 0;
502        if (icasType < 2) ptbinNarrowY = 10;   // 6;
503        else              ptbinNarrowY =  3;   // 2;
504        // - SOME DEFINITIONS
505        Float_t lowlimmass;
506        Float_t uplimmass;
507        Float_t lowgausslim;
508        Float_t upgausslim;
509        if (icasType==0||icasType==1) {
510            lowlimmass=1.30;
511            uplimmass=1.34;
512            lowgausslim=1.312;
513            upgausslim=1.332;
514        } else {
515            lowlimmass=1.645;
516            uplimmass=1.70;
517            lowgausslim=1.668;
518            upgausslim=1.678;
519        }
520        TF1*  fitinvmass = new TF1("fitinvmass","gaus(0)+pol2(3)",lowlimmass,uplimmass);
521        fitinvmass->SetParName(0, "cnstntG");
522        fitinvmass->SetParName(1, "meanG");
523        fitinvmass->SetParName(2, "sigmaG");
524        fitinvmass->SetParLimits(0,0.,500000.);
525        if (icasType==0||icasType==1) {
526            fitinvmass->SetParameter(1, 1.32171);
527            fitinvmass->SetParLimits(1, 1.31,1.33);
528            fitinvmass->SetParLimits(2,0.001,0.005);
529        } else {
530            fitinvmass->SetParameter(1, 1.67245);
531            fitinvmass->SetParLimits(1, 1.664,1.68);
532            fitinvmass->SetParLimits(2,0.0008,0.006);
533        }
534        hvar18->Fit("fitinvmass","rimeN");
535        fitinvmass->SetLineColor(kRed);
536        fitinvmass->Draw("same");
537        Float_t meanGauss   = fitinvmass->GetParameter(1);
538        Float_t sigmaGauss  = fitinvmass->GetParameter(2);
539        cout<<"Mean: "<<meanGauss<<endl;
540        cout<<"Sigma: "<<sigmaGauss<<endl;
541      //Pad2: Text
542       c4->cd(2);
543        Float_t refwidth = 0.002;
544       TPaveText *pave1 = new TPaveText(0.05,0.3,0.95,0.5);
545       pave1->SetFillColor(0);
546       pave1->SetTextSize(0.04);
547       pave1->SetTextAlign(12);
548       if (icasType < 2) pave1->AddText("PDG mass: 1.32171 GeV/c^{2}");
549       else              pave1->AddText("PDG mass: 1.67245 GeV/c^{2}");
550       pave1->AddText(Form("#color[1]{Mass form Fit: %.5f #pm %.5f GeV/c^{2}}",meanGauss,sigmaGauss));
551       if (sigmaGauss > refwidth - 0.0003 && sigmaGauss < refwidth + 0.0003) pave1->AddText("#color[3]{OK!! The width is compatible with standard.}");
552       else                                                                  pave1->AddText("#color[2]{NOT OK!! Problem.}");
553       pave1->Draw();
554       cout<<"   "<<refwidth - 0.0003<<"<"<<sigmaGauss<<"<"<<refwidth + 0.0003<<endl;
555     
556      //DEFINE 5st CANVAS AND DRAW PLOTS
557      if (collidingsystem == 0) {
558        TCanvas *c5 = new TCanvas("c5","" );
559          c5->Divide(2,1);
560             //Pad 1: centrality
561             c5->cd(1);
562             TH1D *hvar16 = cf->ShowProjection(19,icasType);
563             hvar16->Draw("histo");
564
565             hvar16->SetName(Form("fig_lf_multistrange_16", hvar16->GetName()));
566             hvar16->Write();
567
568             //Pad 2: track multiplicity
569             c5->cd(2);
570             TH1D *hvar17 = cf->ShowProjection(20,icasType);
571             hvar17->Draw("histo");
572             
573             hvar17->SetName(Form("fig_lf_multistrange_17", hvar17->GetName()));
574             hvar17->Write();
575             
576             
577             c5->SaveAs(Form("fig_lf_MultistrangeQA_5.%s",suffix.Data()));
578      }
579      
580      fout->Close();
581 }
582
583
584
585
586 Bool_t checkUnderTheLimit(TH1D *lHist, Double_t limit) {
587
588       Int_t binlimit = lHist->FindBin(limit);
589       Bool_t checkOk = kTRUE;
590       for (Int_t i = 1; i < binlimit; i++) {
591            Int_t content = 0;
592            content = lHist->GetBinContent(i);
593            if (content != 0) checkOk = kFALSE;
594            //cout<<"Content bin "<<i<<": "<<content<<endl;
595       }
596       return checkOk;
597
598 }
599
600
601 Bool_t checkOverTheLimit(TH1D *lHist, Double_t limit) {
602
603       Int_t binlimit = lHist->FindBin(limit);
604       Int_t lastbin  = lHist->GetNbinsX();
605       Bool_t checkOk = kTRUE;
606       for (Int_t i = binlimit; i < lastbin+1; i++) {
607            Int_t content = 0;
608            content = lHist->GetBinContent(i);
609            if (content != 0) checkOk = kFALSE;
610            //cout<<"Content bin "<<i<<": "<<content<<endl;
611       }
612       return checkOk;
613
614 }
615
616