]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/QATasks/post/PostProcessQAMultistrange.C
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / QATasks / post / PostProcessQAMultistrange.C
1 //////////////////////////////////////////////////
2 //
3 //  This macro was written by Domenico Colella (domenico.colella@cern.ch).
4 //   -- first version    [12 November 2013]
5 //   -- modified version [3 November 2013]: added inv. mass distr. without PID info
6 //
7 //   ------------------------
8 //   ------ Arguments -------
9 //   ------------------------
10 //   --  icasType          =  0) Xi- 1) Xi+ 2) Omega- 3) Omega+
11 //   --  collidingsystem   =  0) PbPb  1) pp  2) pPb
12 //   --  isMC              =  kTRUE if running on MC production 
13 //   --  fileDir           =  "Input file directory"
14 //   --  filein            =  "Input file name"
15 //
16 //
17 //   -------------------------------------
18 //   ------ QATask output content --------
19 //   -------------------------------------
20 //   The output produced by the QATask is a CFContainer with 4 steps and 21 variables.
21 //   The meaning of each variable within the container are listed here:
22 //   --  0   = Max DCA Cascade Daughters                 pp: 2.0     PbPb: 0.3     pPb: 2.0
23 //   --  1   = Min DCA Bach To PV                        pp: 0.01    PbPb: 0.03    pPb: 0.03
24 //   --  2   = Min Cascade Cosine Of PA                  pp: 0.98    PbPb: 0.999   pPb: 0.95
25 //   --  3   = Min Cascade Radius Fid. Vol.              pp: 0.2     PbPb: 0.9     pPb: 0.4
26 //   --  4   = Window Invariant Mass Lambda              pp: 0.008   PbPb: 0.0008  pPb: 0.010
27 //   --  5   = Max DCA V0 Daughters                      pp: 1.5     PbPb: 1.0     pPb: 2.0
28 //   --  6   = Min V0 Cosine Of PA To PV                 pp: pT dep. PbPb: 0.98    pPb: 0.95
29 //   --  7   = Min V0 Radius Fid. Vol.                   pp: 0.2     PbPb: 0.9     pPb: 1.0
30 //   --  8   = Min DCA V0 To PV                          pp: 0.01    PbPb: 0.05    pPb: 0.05
31 //   --  9   = Min DCA Pos To PV                         pp: 0.05    PbPb: 0.1     pPb: 0.02
32 //   --  10  = Min DCA Neg To PV                         pp: 0.05    PbPb: 0.1     pPb: 0.02
33 //   --  11  = Invariant Mass distribution for Xi
34 //   --  12  = Invariant Mass distribution for Omega
35 //   --  13  = Transverse Momentum distribution
36 //   --  14  = Rapidity distribution for Xi
37 //   --  15  = Rapidity distribution for Omega
38 //   --  16  = Proper length distribution for the cascade
39 //   --  17  = Proper length distribution for the V0
40 //   --  18  = Min V0 Cosine Of PA To Xi Vertex         pp: pT dep. PbPb: pT dep.
41 //   --  19  = Centrality
42 //   --  20  = ESD track multiplicity
43 //   The last two variables are empty in the case proton-proton collisions.
44 //   In case of MC production one more CFContainer is produced, containing infos on the 
45 //   generated particles. As the previous container, this one is composed by 4 steps, one
46 //   for each cascade and 7 variables:
47 //    -- 0   = Total momentum
48 //    -- 1   = Transverse momentum
49 //    -- 2   = Rapidity
50 //    -- 3   = Pseudo-rapidity
51 //    -- 4   = Theta angle
52 //    -- 5   = Phi angle
53 //    -- 6   = Centrality
54 //    The previous container is still produced with the informations from the reconstructed
55 //    particles.
56 //
57 //
58 //   -----------------------------------
59 //   ------ Present Macro Checks -------
60 //   -----------------------------------
61 //   Using this macro many checks on the cascade topological reconstruction procedure
62 //   can be performed. In particular, the shape and the limit for the topological 
63 //   variable distributions as well as other kinematical variable distributions. The
64 //   reconstruction of the cascades are performed using two classes AliCascadeVertexer.cxx 
65 //   and AliV0vertexer.cxx contained in /STEER/ESD/ folder in Aliroot.
66 //   In the following are listed the contents of each page of the produced pdf:
67 //   
68 //   -- [Page 1] Distributions for the variables: 
69 //                DCA cascade daughters,  Bachelor IP to PV, 
70 //                Cascade cosine of PA,   Cascade radius of fiducial volume, 
71 //                Invariant mass Lambda,  DCA V0 daughters.
72 //   -- [Page 2] Distributions for the variables:
73 //                V0 cosine of PA to PV,  Min V0 Radius fiducial volume, 
74 //                Min DCA V0 To PV,       Min DCA positive To PV, 
75 //                Min DCA negative To PV, V0 cosine of PA to XiV
76 //   -- [Page 3] Distributions for the variables;
77 //                InvMass,                Transverse momentum,
78 //                Rapidity,               Cascade proper length, 
79 //                V0 proper length.
80 //   -- [Page 4] Check on the invariant mass distribution fit.
81 //   -- [Page 5] Check on the invariant mass distribution, without the PID info, fit. 
82 //   -- [Page 6] Only in case of PbPb or pPb collisions, the event centrality
83 //               distribution.
84 //   -- [Page 7] Only in case of MC production, distributions for the MC generated
85 //               particles, of the variables:
86 //                Total momentum,         Transverse momentum,
87 //                Rapidity,               Pseudo-rapidity,
88 //                Theta angle,            Phi angle,
89 //
90 //
91 //////////////////////////////////////////////////////
92
93
94
95
96 class AliCFContainer;
97
98 void PostProcessQAMultistrange(Int_t   icasType        = 0,                             // 0) Xi- 1) Xi+ 2) Omega- 3) Omega+
99                                Int_t   collidingsystem = 0,                             // 0) PbPb  1) pp 2) pPb
100                                Bool_t  isMC            = kFALSE,                         // kTRUE-->MC and kFALSE-->Exp.
101                                Char_t *fileDir         = ".",                           // Input file directory
102                                Char_t *filein          = "AnalysisResults.root"         // Input file name
103                               ) {
104
105
106      //___________________
107      //DEFINE DRAW OPTIONS
108      gStyle->SetOptStat(1110);
109      gStyle->SetOptStat(kFALSE);
110      gStyle->SetOptTitle(kFALSE);
111      gStyle->SetFrameLineWidth(2.5);
112      gStyle->SetCanvasColor(0);
113      gStyle->SetPadColor(0);
114      gStyle->SetHistLineWidth(2.5);
115      gStyle->SetLabelSize(0.05, "x");
116      gStyle->SetLabelSize(0.05, "y");
117      gStyle->SetTitleSize(0.05, "x");
118      gStyle->SetTitleSize(0.05, "y");
119      gStyle->SetTitleOffset(1.1, "x");
120      gStyle->SetPadBottomMargin(0.14);
121
122      //_______________________
123      //SOURCE USEFUL LIBRARIES
124      gSystem->Load("libANALYSIS");
125      gSystem->Load("libANALYSISalice");
126      gSystem->Load("libCORRFW");
127
128      //_________________________________
129      //SOURCE THE FILE AND THE CONTAINER
130      TFile *f1 = new TFile(Form("%s/%s",fileDir,filein));
131      AliCFContainer *cf = (AliCFContainer*) (f1->Get("PWGLFStrangeness.outputCheckCascade/fCFContCascadeCuts"));  
132      TList *hlist = (TList*) f1->Get("PWGLFStrangeness.outputCheckCascade/fListHistMultistrangeQA");
133  
134      //____________
135      //DEEFINE TEXT
136      TLatex* t1 = new TLatex(0.6,0.55,"#color[3]{OK!!}");
137      t1->SetTextSize(0.1);
138      t1->SetNDC();
139      TLatex* t2 = new TLatex(0.6,0.55,"#color[2]{NOT OK!!}");
140      t2->SetTextSize(0.1);
141      t2->SetNDC();
142      t2->SetTextColor(2);
143      TLatex* tcasc;
144      if      (icasType == 0) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{-}}");
145      else if (icasType == 1) tcasc = new TLatex(0.8,0.7,"#color[1]{#Xi^{+}}");
146      else if (icasType == 2) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{-}}");
147      else if (icasType == 3) tcasc = new TLatex(0.8,0.7,"#color[1]{#Omega^{+}}");
148      tcasc->SetTextSize(0.15);
149      tcasc->SetNDC();
150      tcasc->SetTextColor(2);
151      TLatex* tpdgmass;
152      if      (icasType == 0) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
153      else if (icasType == 1) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.321 GeV/c^{2}}");
154      else if (icasType == 2) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
155      else if (icasType == 3) tpdgmass = new TLatex(0.55,0.7,"#color[1]{PDG mass: 1.672 GeV/c^{2}}");
156      tpdgmass->SetTextSize(0.07);
157      tpdgmass->SetNDC();
158      tpdgmass->SetTextColor(2);
159
160      
161      //________________________________ 
162      //DEFINE 1st CANVAS AND DRAW PLOTS
163      TCanvas *c1 = new TCanvas("c1","",1200,800);
164      c1->Divide(2,3); 
165        //Pad 1: DCA cascade daughters
166        c1->cd(1);
167        gPad->SetLogy();
168        TH1D *hvar0 = cf->ShowProjection(0,icasType);
169        hvar0->Draw("histo");
170        Double_t x0;
171        if      (collidingsystem == 0) x0 = 0.3;
172        else if (collidingsystem == 1) x0 = 2.0;
173        else if (collidingsystem == 2) x0 = 2.0;
174        TLine *line0 = new TLine(x0,0.,x0,hvar0->GetBinContent(hvar0->GetMaximumBin()));
175        line0->SetLineColor(kRed);
176        line0->SetLineStyle(9);
177        line0->SetLineWidth(2.0);
178        line0->Draw("same");
179           Bool_t check_0 = checkOverTheLimit(hvar0, x0);
180           if (check_0) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
181           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
182        tcasc->Draw();
183        //Pad 2: Bachelor IP to PV
184        c1->cd(2);
185        gPad->SetLogy();
186        TH1D *hvar1 = cf->ShowProjection(1,icasType);
187        hvar1->GetXaxis()->SetRangeUser(0.,0.24);
188        hvar1->Draw("histo");
189        Double_t x1;
190        if      (collidingsystem == 0) x1 = 0.03;
191        else if (collidingsystem == 1) x1 = 0.01;
192        else if (collidingsystem == 2) x1 = 0.03;
193        TLine *line1 = new TLine(x1,0.,x1,hvar1->GetBinContent(hvar1->GetMaximumBin()));
194        line1->SetLineColor(kRed);
195        line1->SetLineStyle(9);
196        line1->SetLineWidth(2.0);
197        line1->Draw("same");
198           Bool_t check_1 = checkUnderTheLimit(hvar1, x1);
199           if (check_1) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
200           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
201        //Pad 3: Cascade cosine of Pointing Angle
202        c1->cd(3);
203        gPad->SetLogy();
204        TH1D *hvar2 = cf->ShowProjection(2,icasType);
205        Double_t max2 = hvar2->GetBinContent(hvar2->GetMaximumBin());
206        hvar2->GetYaxis()->SetRangeUser(0.01,max2*1.5);
207        hvar2->Draw("histo");
208        Double_t x2;
209        if      (collidingsystem == 0) x2 = 0.999;
210        else if (collidingsystem == 1) x2 = 0.98;
211        else if (collidingsystem == 2) x2 = 0.95;
212        TLine *line2 = new TLine(x2,0.,x2,hvar2->GetBinContent(hvar2->GetMaximumBin()));
213        line2->SetLineColor(kRed);
214        line2->SetLineStyle(9);
215        line2->SetLineWidth(2.0);
216        line2->Draw("same");
217        line1->Draw("same");
218           Bool_t check_2 = checkUnderTheLimit(hvar2, x2);
219           if (check_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
220           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
221        //Pad 4: Cascade radius of fiducial volume
222        c1->cd(4);
223        gPad->SetLogy();
224        TH1D *hvar3 = cf->ShowProjection(3,icasType);
225        hvar3->GetXaxis()->SetRangeUser(0.,3.8);
226        hvar3->Draw("histo");
227        Double_t x3;
228        if      (collidingsystem == 0) x3 = 0.9;
229        else if (collidingsystem == 1) x3 = 0.2;
230        else if (collidingsystem == 2) x3 = 0.4;
231        TLine *line3 = new TLine(x3,0.,x3,hvar3->GetBinContent(hvar3->GetMaximumBin()));
232        line3->SetLineColor(kRed);
233        line3->SetLineStyle(9);
234        line3->SetLineWidth(2.0);
235        line3->Draw("same");
236           Bool_t check_3 = checkUnderTheLimit(hvar3, x3);
237           if (check_3) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
238           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
239        //Pad 5: Invariant mass Lambda
240        c1->cd(5);
241        TH1D *hvar4 = cf->ShowProjection(4,icasType);
242        hvar4->Draw("histo");
243        Double_t x41;
244        if      (collidingsystem < 2)  x41 = 1.116 + 0.008;
245        else if (collidingsystem == 2) x41 = 1.116 + 0.010;
246        TLine *line41 = new TLine(x41,0.,x41,hvar4->GetBinContent(hvar4->GetMaximumBin()));
247        line41->SetLineColor(kRed);
248        line41->SetLineStyle(9);
249        line41->SetLineWidth(2.0);
250        line41->Draw("same");
251        Double_t x42;
252        if      (collidingsystem < 2)  x42 = 1.115 - 0.008;
253        else if (collidingsystem == 2) x42 = 1.115 - 0.010;
254        TLine *line42 = new TLine(x42,0.,x42,hvar4->GetBinContent(hvar4->GetMaximumBin()));
255        line42->SetLineColor(kRed);
256        line42->SetLineStyle(9);
257        line42->SetLineWidth(2.0);
258        line42->Draw("same");
259           Bool_t check_4_1 = checkUnderTheLimit(hvar4, x42);
260           Bool_t check_4_2 = checkOverTheLimit(hvar4, x41);
261           if (check_4_1 && check_4_2) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
262           else                        { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
263        //Pad 6: DCA V0 daughters
264        c1->cd(6);
265        gPad->SetLogy();
266        TH1D *hvar5 = cf->ShowProjection(5,icasType);
267        hvar5->Draw("histo");
268        Double_t x5;
269        if      (collidingsystem == 0) x5 = 1.0;
270        else if (collidingsystem == 1) x5 = 1.5;
271        else if (collidingsystem == 2) x5 = 2.0;
272        TLine *line5 = new TLine(x5,0.,x5,hvar5->GetBinContent(hvar5->GetMaximumBin()));
273        line5->SetLineColor(kRed);
274        line5->SetLineStyle(9);
275        line5->SetLineWidth(2.0);
276        line5->Draw("same");
277           Bool_t check_5 = checkOverTheLimit(hvar5, x5);
278           if (check_5) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
279           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
280      c1->SaveAs("fig_lf_Multistrange.pdf(");
281     
282      //________________________________
283      //DEFINE 2nd CANVAS AND DRAW PLOTS
284      TCanvas *c2 = new TCanvas("c2","",1200,800);
285      c2->Divide(2,3);
286        //Pad 1: V0 cosine of Pointing Angle to PV
287        c2->cd(1);
288        gPad->SetLogy();
289        TH1D *hvar6 = cf->ShowProjection(6,icasType);
290        Double_t max6 = hvar6->GetBinContent(hvar6->GetMaximumBin());
291        hvar6->GetYaxis()->SetRangeUser(0.01,max6*1.5);
292        hvar6->Draw("histo");
293        tcasc->Draw();
294        //Pad 2: Min V0 Radius Fid. Vol.  
295        c2->cd(2);
296        gPad->SetLogy();
297        TH1D *hvar7 = cf->ShowProjection(7,icasType);
298        hvar7->GetXaxis()->SetRangeUser(0.,3.0);
299        hvar7->Draw("histo");
300        Double_t x7;
301        if      (collidingsystem == 0) x7 = 0.9;
302        else if (collidingsystem == 1) x7 = 0.2;
303        else if (collidingsystem == 2) x7 = 0.4;
304        TLine *line7 = new TLine(x7,0.,x7,hvar7->GetBinContent(hvar7->GetMaximumBin()));
305        line7->SetLineColor(kRed);
306        line7->SetLineStyle(9);
307        line7->SetLineWidth(2.0);
308        line7->Draw("same");
309           Bool_t check_7 = checkUnderTheLimit(hvar7, x7);
310           if (check_7) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
311           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
312        //Pad3: Min DCA V0 To PV
313        c2->cd(3);
314        gPad->SetLogy();
315        TH1D *hvar8 = cf->ShowProjection(8,icasType);
316        hvar8->GetXaxis()->SetRangeUser(0.,0.3);
317        hvar8->Draw("histo");
318        Double_t x8;
319        if      (collidingsystem == 0) x8 = 0.05;
320        else if (collidingsystem == 1) x8 = 0.01;
321        else if (collidingsystem == 2) x8 = 0.05;
322        TLine *line8 = new TLine(x8,0.,x8,hvar8->GetBinContent(hvar8->GetMaximumBin()));
323        line8->SetLineColor(kRed);
324        line8->SetLineStyle(9);
325        line8->SetLineWidth(2.0);
326        line8->Draw("same");
327           Bool_t check_8 = checkUnderTheLimit(hvar8, x8);
328           if (check_8) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
329           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
330        //Pad 4: Min DCA Pos To PV
331        c2->cd(4);
332        gPad->SetLogy();
333        TH1D *hvar9 = cf->ShowProjection(9,icasType);
334        hvar9->GetXaxis()->SetRangeUser(0.,0.2);
335        hvar9->Draw("histo");
336        Double_t x9;
337        if      (collidingsystem == 0) x9 = 0.1;
338        else if (collidingsystem == 1) x9 = 0.05;
339        else if (collidingsystem == 2) x9 = 0.02;
340        TLine *line9 = new TLine(x9,0.,x9,hvar9->GetBinContent(hvar9->GetMaximumBin()));
341        line9->SetLineColor(kRed);
342        line9->SetLineStyle(9);
343        line9->SetLineWidth(2.0);
344        line9->Draw("same");
345           Bool_t check_9 = checkUnderTheLimit(hvar9, x9);
346           if (check_9) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
347           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
348        //Pad 5: Min DCA Neg To PV
349        c2->cd(5);
350        gPad->SetLogy();
351        TH1D *hvar10 = cf->ShowProjection(10,icasType);
352        hvar10->GetXaxis()->SetRangeUser(0.,0.2);
353        hvar10->Draw("histo");
354        Double_t x10;
355        if      (collidingsystem == 0) x10 = 0.1;
356        else if (collidingsystem == 1) x10 = 0.05;
357        else if (collidingsystem == 2) x10 = 0.02;
358        TLine *line10 = new TLine(x10,0.,x10,hvar10->GetBinContent(hvar10->GetMaximumBin()));
359        line10->SetLineColor(kRed);
360        line10->SetLineStyle(9);
361        line10->SetLineWidth(2.0);
362        line10->Draw("same");
363           Bool_t check_10 = checkUnderTheLimit(hvar10, x10);
364           if (check_10) { cout<<"The cut is OK!!"<<endl; t1->Draw(); }
365           else         { cout<<"The cut is NOT OK!!"<<endl; t2->Draw(); }
366        //Pad 6: V0 cosine of Pointing Angle to Xi vtx
367        c2->cd(6);
368        gPad->SetLogy();
369        TH1D *hvar20 = cf->ShowProjection(18,icasType);
370        Double_t max20 = hvar20->GetBinContent(hvar20->GetMaximumBin());
371        hvar20->GetYaxis()->SetRangeUser(0.01,max20*1.5);
372        hvar20->Draw("histo");
373      c2->SaveAs("fig_lf_Multistrange.pdf");
374
375      //________________________________
376      //DEFINE 3rd CANVAS AND DRAW PLOTS
377      TCanvas *c3 = new TCanvas("c3","",1200,800);
378      c3->Divide(2,2);
379        //Pad 1: Transverse momentum
380        c3->cd(1);
381        TH1D *hvar13 = cf->ShowProjection(13,icasType);
382        hvar13->Draw("histo");
383        tcasc->Draw();
384        //Pad 2: Y
385        c3->cd(2);
386        TH1D *hvar14 = cf->ShowProjection(14+icasType/2,icasType);
387        hvar14->Draw("histo");
388        //Pad 3: Cascade proper length
389        c3->cd(3);
390        TH1D *hvar18;
391        hvar18 = cf->ShowProjection(16,icasType);
392        hvar18->GetXaxis()->SetRangeUser(0.,90.);
393        hvar18->Draw("histo");
394        //Pad 4: V0 proper length 
395        c3->cd(4);
396        TH1D *hvar19;
397        hvar19 = cf->ShowProjection(17,icasType);
398        hvar19->GetXaxis()->SetRangeUser(0.,90.);
399        hvar19->Draw("histo");
400        //Pad 6
401        // empty 
402      c3->SaveAs("fig_lf_Multistrange.pdf");
403
404      //________________________________ 
405      //DEFINE 4th CANVAS AND DRAW PLOTS
406      TCanvas *c4 = new TCanvas("c4","",1200,800);
407      c4->Divide(2,2);
408        //Pad1: invariant mass fit
409        c4->cd(1);
410        TH1D *hvar23 = cf->ShowProjection(11+icasType/2,icasType);
411        hvar23->Draw("histo");
412        tcasc->Draw();
413         // - SOME PARAMETER VALUE
414         Bool_t kfitgauss = kFALSE;
415         Bool_t kfitleft  = kFALSE;
416         Bool_t kfitright = kFALSE;
417         // - SOME DEFINITIONS
418         Float_t lowlimmass;
419         Float_t uplimmass;
420         Float_t lowgausslim;
421         Float_t upgausslim;
422         if (icasType==0||icasType==1) {
423             lowlimmass=1.30;
424             uplimmass=1.34;
425             lowgausslim=1.312;
426             upgausslim=1.332;
427         } else {
428             lowlimmass=1.645;
429             uplimmass=1.70;
430             lowgausslim=1.668;
431             upgausslim=1.678;
432         }
433         TF1*  fitinvmass = new TF1("fitinvmass","gaus(0)+pol2(3)",lowlimmass,uplimmass);
434         fitinvmass->SetParName(0, "cnstntG");
435         fitinvmass->SetParName(1, "meanG");
436         fitinvmass->SetParName(2, "sigmaG");
437         fitinvmass->SetParLimits(0,0.,500000.);
438         if (icasType==0||icasType==1) {
439             fitinvmass->SetParameter(1, 1.32171);
440             fitinvmass->SetParLimits(1, 1.31,1.33);
441             fitinvmass->SetParLimits(2,0.001,0.005);
442         } else {
443             fitinvmass->SetParameter(1, 1.67245);
444             fitinvmass->SetParLimits(1, 1.664,1.68);
445             fitinvmass->SetParLimits(2,0.0008,0.006);
446         }
447         hvar23->Fit("fitinvmass","rimeN");
448         fitinvmass->SetLineColor(kRed);
449         fitinvmass->Draw("same");
450         Float_t meanGauss   = fitinvmass->GetParameter(1);
451         Float_t sigmaGauss  = fitinvmass->GetParameter(2);
452        cout<<"Mean: "<<meanGauss<<endl;
453        cout<<"Sigma: "<<sigmaGauss<<endl;
454        //Pad2: Text
455        c4->cd(2);
456        Float_t refwidth = 0.002;
457        if (icasType > 1) refwidth = 0.0025;
458        TPaveText *pave1 = new TPaveText(0.05,0.3,0.95,0.7);
459        pave1->SetFillColor(0);
460        pave1->SetTextSize(0.04);
461        pave1->SetTextAlign(12);
462        pave1->AddText("Invariant mass distribution #color[4]{WITH} PID on dauther tracks");
463        if (icasType < 2) pave1->AddText("PDG mass: 1.32171 GeV/c^{2}");
464        else              pave1->AddText("PDG mass: 1.67245 GeV/c^{2}");
465        pave1->AddText(Form("#color[1]{Mass from Fit: %.5f #pm %.5f GeV/c^{2}}",meanGauss,sigmaGauss));
466        if (sigmaGauss > refwidth - 0.0003 && sigmaGauss < refwidth + 0.0003) pave1->AddText("#color[3]{OK!! The width is compatible with standard.}");
467        else                                                                  pave1->AddText("#color[2]{NOT OK!! Problem.}");
468        pave1->Draw();
469        cout<<"   "<<refwidth - 0.0003<<"<"<<sigmaGauss<<"<"<<refwidth + 0.0003<<endl;
470        //Pad3: invariant mass fit
471        c4->cd(3);
472        TH1D *hvar24;
473        if      (icasType == 0) hvar24 = (TH1D*) hlist->FindObject(Form("fHistMassXiMinus"));
474        else if (icasType == 1) hvar24 = (TH1D*) hlist->FindObject(Form("fHistMassXiPlus"));
475        else if (icasType == 2) hvar24 = (TH1D*) hlist->FindObject(Form("fHistMassOmegaMinus"));
476        else if (icasType == 3) hvar24 = (TH1D*) hlist->FindObject(Form("fHistMassOmegaPlus"));
477        hvar24->Draw("histo");
478         TF1*  fitinvmassnoPID = new TF1("fitinvmassnoPID","gaus(0)+pol2(3)",lowlimmass,uplimmass);
479         fitinvmassnoPID->SetParName(0, "cnstntG");
480         fitinvmassnoPID->SetParName(1, "meanG");
481         fitinvmassnoPID->SetParName(2, "sigmaG");
482         fitinvmassnoPID->SetParLimits(0,0.,500000.);
483         if (icasType==0||icasType==1) {
484             fitinvmassnoPID->SetParameter(1, 1.32171);
485             fitinvmassnoPID->SetParLimits(1, 1.31,1.33);
486             fitinvmassnoPID->SetParLimits(2,0.001,0.005);
487         } else {
488             fitinvmassnoPID->SetParameter(1, 1.67245);
489             fitinvmassnoPID->SetParLimits(1, 1.664,1.68);
490             fitinvmassnoPID->SetParLimits(2,0.0008,0.006);
491         }
492         hvar24->Fit("fitinvmassnoPID","rimeN");
493         fitinvmassnoPID->SetLineColor(kRed);
494         fitinvmassnoPID->Draw("same");
495         Float_t meanGaussnoPID   = fitinvmassnoPID->GetParameter(1);
496         Float_t sigmaGaussnoPID  = fitinvmassnoPID->GetParameter(2);
497        cout<<"Mean: "<<meanGaussnoPID<<endl;
498        cout<<"Sigma: "<<sigmaGaussnoPID<<endl;
499        //Pad4: Text
500        c4->cd(4);
501        Float_t refwidth = 0.002;
502        if (icasType > 1) refwidth = 0.0025;
503        TPaveText *pave2 = new TPaveText(0.05,0.3,0.95,0.7);
504        pave2->SetFillColor(0);
505        pave2->SetTextSize(0.04);
506        pave2->SetTextAlign(12);
507        pave2->AddText("Invariant mass distribution #color[4]{WITHOUT} PID on dauther tracks");
508        if (icasType < 2) pave2->AddText("PDG mass: 1.32171 GeV/c^{2}");
509        else              pave2->AddText("PDG mass: 1.67245 GeV/c^{2}");
510        pave2->AddText(Form("#color[1]{Mass from Fit: %.5f #pm %.5f GeV/c^{2}}",meanGaussnoPID,sigmaGaussnoPID));
511        if (sigmaGaussnoPID > refwidth - 0.0003 && sigmaGaussnoPID < refwidth + 0.0003) pave2->AddText("#color[3]{OK!! The width is compatible with standard.}");
512        else                                                                            pave2->AddText("#color[2]{NOT OK!! Problem.}");
513        pave2->Draw();
514        cout<<"   "<<refwidth - 0.0003<<"<"<<sigmaGaussnoPID<<"<"<<refwidth + 0.0003<<endl;
515      c4->SaveAs("fig_lf_Multistrange.pdf");
516
517
518      //________________________________ 
519      //DEFINE 6th CANVAS AND DRAW PLOTS
520      TCanvas *c6 = new TCanvas("c6","",1200,400);//1200,270);
521      c6->Divide(2,1);
522        //Pad 1: Event Selection
523        c6->cd(1);
524        TH1D* hvar26 = hlist->FindObject(Form("fHistEventSel"));
525        hvar26->Draw("histo");
526        //Pad 2: Centrality/Multiplicity distribution
527        c6->cd(2);
528        TH1D *hvar16 = cf->ShowProjection(19,icasType);
529        hvar16->Draw("histo");
530      if      (!isMC) c6->SaveAs("fig_lf_Multistrange.pdf)");
531      else if (isMC) c6->SaveAs("fig_lf_Multistrange.pdf");
532
533
534      //_______________________________
535      //CHECK ON MONTE CARLO PRODUCTION
536      if (isMC) { 
537            
538             AliCFContainer *cfMC = (AliCFContainer*) (f1->Get("PWGLFStrangeness.outputCheckCascade/fCFContCascadeMCgen"));
539             //DEFINE 6th CANVAS AND DRAW PLOTS
540             TCanvas *c7 = new TCanvas("c7","",1200,800);
541             c7->Divide(2,3);
542             //Pad 1: Total Momentum
543             c7->cd(1);
544             TH1D *hvar17 = cfMC->ShowProjection(0,icasType);
545             hvar17->Draw("histo");
546             tcasc->Draw();
547             //Pad 2: Transverse Momentum
548             c7->cd(2);
549             TH1D *hvar18 = cfMC->ShowProjection(1,icasType);
550             hvar18->Draw("histo");
551             //Pad 3: Rapidity (y)
552             c7->cd(3);
553             TH1D *hvar19 = cfMC->ShowProjection(2,icasType);
554             hvar19->Draw("histo");
555             //Pad 4: Pseudo-rapidity (eta)
556             c7->cd(4);
557             TH1D *hvar20 = cfMC->ShowProjection(3,icasType);
558             hvar20->Draw("histo");
559             //Pad 5: Theta
560             c7->cd(5);
561             TH1D *hvar21 = cfMC->ShowProjection(4,icasType);
562             hvar21->Draw("histo");
563             //Pad 6: Phi
564             c7->cd(6);
565             TH1D *hvar22 = cfMC->ShowProjection(5,icasType);
566             hvar22->Draw("histo");
567             c7->SaveAs("fig_lf_Multistrange.pdf)");
568      }
569
570
571 }
572
573
574
575 //______________________
576 Bool_t checkUnderTheLimit(TH1D *lHist, Double_t limit) {
577
578       Int_t binlimit = lHist->FindBin(limit);
579       Bool_t checkOk = kTRUE;
580       for (Int_t i = 1; i < binlimit; i++) {
581            Int_t content = 0;
582            content = lHist->GetBinContent(i);
583            if (content != 0) checkOk = kFALSE;
584            //cout<<"Content bin "<<i<<": "<<content<<endl;
585       }
586       return checkOk;
587
588 }
589
590 //______________________
591 Bool_t checkOverTheLimit(TH1D *lHist, Double_t limit) {
592
593       Int_t binlimit = lHist->FindBin(limit);
594       Int_t lastbin  = lHist->GetNbinsX();
595       Bool_t checkOk = kTRUE;
596       for (Int_t i = binlimit; i < lastbin+1; i++) {
597            Int_t content = 0;
598            content = lHist->GetBinContent(i);
599            if (content != 0) checkOk = kFALSE;
600            //cout<<"Content bin "<<i<<": "<<content<<endl;
601       }
602       return checkOk;
603
604 }
605
606
607