]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/Nuclei/masses/AliAnalysisNucleiMass.cxx
0a4a85ec8789e1c8d70453d3965c3c39abda93b1
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / Nuclei / masses / AliAnalysisNucleiMass.cxx
1
2 #include "AliAnalysisNucleiMass.h"
3
4 // ROOT includes
5 #include <TMath.h>
6 #include "TChain.h"
7
8 // AliRoot includes
9 #include "AliInputEventHandler.h"
10 #include "AliAODEvent.h"
11 #include "AliESDEvent.h"
12 #include "AliVEvent.h"
13 #include "AliAODTrack.h"
14 #include "AliAODPid.h"
15 #include "AliCentrality.h"
16 #include "TH2F.h"
17 #include "TH2D.h"
18 #include "TH1F.h"
19 #include "TF1.h"
20 #include "TGraph.h"
21 #include "TProfile.h"
22 #include "AliESDtrackCuts.h"
23 #include "AliAnalysisManager.h"
24
25 ClassImp(AliAnalysisNucleiMass)
26
27 //_____________________________________________________________________________
28 AliAnalysisNucleiMass::AliAnalysisNucleiMass():
29   AliAnalysisTaskSE(),
30   fMC(kFALSE),
31   FilterBit(16),
32   NminTPCcluster(0),
33   DCAzCUT(100),
34   DCAxyCUT(0.1),
35   kTPCcut(kTRUE),
36   kTPC(0),
37   kTOF(0),
38   fAOD(NULL),
39   fESD(NULL),
40   fEvent(NULL),
41   fPIDResponse(NULL),
42   fList1(new TList()),
43   hNevent(NULL),
44   hZvertex(NULL),    
45   fBetaTofVSp(NULL),
46   hTOFSignalPion(NULL)
47 {
48    // Default constructor (should not be used)
49   fList1->SetName("results");
50 }
51 //______________________________________________________________________________
52 AliAnalysisNucleiMass::AliAnalysisNucleiMass(const char *name):
53   AliAnalysisTaskSE(name),
54   fMC(kFALSE),
55   FilterBit(16),
56   NminTPCcluster(0),
57   DCAzCUT(100),
58   DCAxyCUT(0.1),
59   kTPCcut(kTRUE),
60   kTPC(0),
61   kTOF(0),
62   fAOD(NULL), 
63   fESD(NULL),
64   fEvent(NULL),
65   fPIDResponse(NULL),
66   fList1(new TList()),
67   hNevent(NULL),
68   hZvertex(NULL),    
69   fBetaTofVSp(NULL),
70   hTOFSignalPion(NULL)
71 {
72   DefineOutput(1, TList::Class());
73   fList1->SetName("results");
74 }
75 //_____________________________________________________________________________
76 AliAnalysisNucleiMass::~AliAnalysisNucleiMass()
77 {
78   if(fList1) delete fList1;
79 }
80 //______________________________________________________________________________
81 void AliAnalysisNucleiMass::UserCreateOutputObjects()
82 {
83   
84   hNevent = new TH1F("hNevent","Centrality",20,0,100);
85
86   hZvertex = new TH1F("hZvertex","Vertex distribution of accepted events; z vertex (cm)",120,-15,15);
87
88   hTOFSignalPion = new TH1F("hTOFSignalPion","TOF signal 0.9<p_{T}<1.0; t-t_{exp}^{#pi} (ps)",1500,-1500,1500);
89   
90   fdEdxVSp[0] = new TH2F("fdEdxVSp","dE/dx vs p; p/|Z| (GeV/c); dE/dx_{TPC} (a.u.)",500,0,5,2000,0,1000);
91   fdEdxVSp[1] = new TH2F("fdEdxVSp_pos","dE/dx vs p positive charge; p/|Z| (GeV/c); dE/dx_{TPC} (a.u.)",500,0,5,2000,0,1000);
92   fdEdxVSp[2] = new TH2F("fdEdxVSp_neg","dE/dx vs p negative charge; p/|Z| (GeV/c); dE/dx_{TPC} (a.u.)",500,0,5,2000,0,1000);
93   
94   fBetaTofVSp = new TH2F("fBetaTofVSp","#beta_{TOF} vs p; p(GeV/c); #beta_{TOF}",1000,0,5,1300,0.4,1.05);
95   
96   fM2vsP_NoTpcCut[0] = new TH2F("fM2vsP_NoTpcCut","M_{TOF}^{2} vs p; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
97   fM2vsP_NoTpcCut[1] = new TH2F("fM2vsP_NoTpcCut_Positive","M_{TOF}^{2} vs p Pos Part; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
98   fM2vsP_NoTpcCut[2] = new TH2F("fM2vsP_NoTpcCut_Negative","M_{TOF}^{2} vs p Neg Part; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
99   
100   fM2vsP_NoTpcCut_DCAxyCut[0] = new TH2F("fM2vsP_NoTpcCut_DCAxycut","M_{TOF}^{2} vs p with DCAxy cut; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
101   fM2vsP_NoTpcCut_DCAxyCut[1] = new TH2F("fM2vsP_NoTpcCut_Positive_DCAxycut","M_{TOF}^{2} vs p Pos Part with DCAxy cut; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
102   fM2vsP_NoTpcCut_DCAxyCut[2] = new TH2F("fM2vsP_NoTpcCut_Negative_DCAxycut","M_{TOF}^{2} vs p Neg Part with DCAxy cut; M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4}); p/|Z| (GeV/c)",8000,0,10,200,0,10);
103
104   fM2vsZ[0] = new TH2F("fM2vsZ","M_{TOF}^{2} vs Z^{2} Integrated p_{T};Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
105   fM2vsZ[1] = new TH2F("fM2vsZ_0.3pT0.5","M_{TOF}^{2} vs Z^{2} 0.3<pT<0.5;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
106   fM2vsZ[2] = new TH2F("fM2vsZ_0.5pT1.0","M_{TOF}^{2} vs Z^{2} 0.5<pT<1.0;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
107   fM2vsZ[3] = new TH2F("fM2vsZ_1.0pT1.5","M_{TOF}^{2} vs Z^{2} 1.0<pT<1.5;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
108   fM2vsZ[4] = new TH2F("fM2vsZ_1.5pT2.0","M_{TOF}^{2} vs Z^{2} 1.5<pT<2.0;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
109   fM2vsZ[5] = new TH2F("fM2vsZ_2.0pT2.5","M_{TOF}^{2} vs Z^{2} 2.0<pT<2.5;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
110   fM2vsZ[6] = new TH2F("fM2vsZ_2.5pT3.0","M_{TOF}^{2} vs Z^{2} 2.5<pT<3.0;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
111   fM2vsZ[7] = new TH2F("fM2vsZ_3.0pT3.5","M_{TOF}^{2} vs Z^{2} 3.0<pT<3.5;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
112   fM2vsZ[8] = new TH2F("fM2vsZ_3.5pT4.0","M_{TOF}^{2} vs Z^{2} 3.5<pT<4.0;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
113   fM2vsZ[9] = new TH2F("fM2vsZ_4.0pT4.5","M_{TOF}^{2} vs Z^{2} 4.0<pT<4.5;Z;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4})",4000,-4,4,1000,0,10);
114   
115   
116   char namePart[9][30];
117   char namePart_par_TPC[9][40];
118   char namePart_title_TPC[9][120];
119   
120   char namePart_par_TOF[9][40];
121   char namePart_title_TOF[9][120];
122
123   char namePart_par_ProfileTPC[9][40];
124   char namePart_title_ProfileTPC[9][80];
125
126   char namePart_par_ProfileTOF[9][40];
127   char namePart_title_ProfileTOF[9][80];
128  
129   
130   snprintf(namePart[0],20,"e");
131   snprintf(namePart[1],20,"#mu");
132   snprintf(namePart[2],20,"#pi");
133   snprintf(namePart[3],20,"K");
134   snprintf(namePart[4],20,"p");
135   snprintf(namePart[5],20,"d");
136   snprintf(namePart[6],20,"t");
137   snprintf(namePart[7],20,"He3");
138   snprintf(namePart[8],20,"He4");
139   
140   for(Int_t i=0;i<9;i++) {
141     snprintf(namePart_par_TPC[i],40,"NsigmaTPC_%s",namePart[i]);
142     snprintf(namePart_title_TPC[i],120,"NsigmaTPC_%s;p_{T} (GeV/c);n_{#sigma_{TPC}}^{%s}",namePart[i],namePart[i]);
143   
144     snprintf(namePart_par_TOF[i],40,"NsigmaTOF_%s",namePart[i]);
145     snprintf(namePart_title_TOF[i],120,"NsigmaTOF_%s;p_{T} (GeV/c);n_{#sigma_{TOF}}^{%s}",namePart[i],namePart[i]);
146  
147     snprintf(namePart_par_ProfileTPC[i],40,"hDeDxExp_%s",namePart[i]);
148     snprintf(namePart_title_ProfileTPC[i],80,"hDeDxExp_%s;p (GeV/c);dE/dx_{TPC} (a.u.)",namePart[i]);
149  
150     snprintf(namePart_par_ProfileTOF[i],40,"hBetaVsP_Exp_%s",namePart[i]);
151     snprintf(namePart_title_ProfileTOF[i],80,"hBetaVsP_Exp%s;p (GeV/c); #beta_{TOF}",namePart[i]);
152   }
153   
154   for(Int_t i=0;i<9;i++) {
155     fNsigmaTPC[i] = new TH2F(namePart_par_TPC[i],namePart_title_TPC[i],250,0,5,200,-5,5);
156     fNsigmaTPC[i]->GetYaxis()->CenterTitle();
157     fNsigmaTOF[i] = new TH2F(namePart_par_TOF[i],namePart_title_TOF[i],250,0,5,200,-5,5);
158     fNsigmaTOF[i]->GetYaxis()->CenterTitle();
159     hDeDxExp[i] = new TProfile(namePart_par_ProfileTPC[i],namePart_title_ProfileTPC[i],500,0,5,0,1000,"");
160     hBetaExp[i] = new TProfile(namePart_par_ProfileTOF[i],namePart_title_ProfileTOF[i],400,0,5,0.4,1.05,"");
161   }
162   
163   char namePart_par_TPCvsP_kTOFtrue[18][80];
164   char namePart_title_TPCvsP_kTOFtrue[18][120];
165
166   char name[18][30];
167
168   snprintf(name[0],20,"e^{+}");
169   snprintf(name[1],20,"#mu^{+}");
170   snprintf(name[2],20,"#pi^{+}");
171   snprintf(name[3],20,"K^{+}");
172   snprintf(name[4],20,"p");
173   snprintf(name[5],20,"d");
174   snprintf(name[6],20,"t");
175   snprintf(name[7],20,"He3");
176   snprintf(name[8],20,"He4");
177   
178   snprintf(name[9],20,"e^{-}");
179   snprintf(name[10],20,"#mu^{-}");
180   snprintf(name[11],20,"#pi^{-}");
181   snprintf(name[12],20,"K^{-}");
182   snprintf(name[13],20,"#bar{p}");
183   snprintf(name[14],20,"#bar{d}");
184   snprintf(name[15],20,"#bar{t}");
185   snprintf(name[16],20,"#bar{He3}");
186   snprintf(name[17],20,"#bar{He4}");
187   
188   for(Int_t iS=0;iS<18;iS++) {
189     snprintf(namePart_par_TPCvsP_kTOFtrue[iS],40,"NsigmaTPCvsP_kTOFout&&kTIME_%s",name[iS]);
190     snprintf(namePart_title_TPCvsP_kTOFtrue[iS],150,"NsigmaTPCvsP_kTOFout&&kTIME_%s;p (GeV/c);n_{#sigma_{TPC}}^{%s}",name[iS],name[iS]);
191   }
192   
193   for(Int_t iS=0;iS<18;iS++) {
194     fNsigmaTPCvsP_kTOFtrue[iS] = new TH2F(namePart_par_TPCvsP_kTOFtrue[iS],namePart_title_TPCvsP_kTOFtrue[iS],250,0,5,200,-5,5);
195     fNsigmaTPCvsP_kTOFtrue[iS]->GetYaxis()->CenterTitle();
196   }
197
198   char name_par_MvsP[18][60];
199   char name_title_MvsP[18][150];
200   
201   for(Int_t i=0;i<18;i++) {
202     snprintf(name_par_MvsP[i],60,"M2vsP_%s",name[i]);
203     snprintf(name_title_MvsP[i],150,"M_{TOF}^{2}_%s_2#sigma_TPCcut;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4});p/|Z| (GeV/c)",name[i]);
204   }
205
206   for (Int_t i=0;i<18;i++) fM2vsP[i] = new TH2F(name_par_MvsP[i],name_title_MvsP[i],8000,0,10,200,0,10);
207   
208   char name_par_MvsP_DCAxyCut[18][60];
209   char name_title_MvsP_DCAxyCut[18][150];
210   
211   for(Int_t i=0;i<18;i++) {
212     snprintf(name_par_MvsP_DCAxyCut[i],60,"M2vsP_DCAxyCut_%s",name[i]);
213     snprintf(name_title_MvsP_DCAxyCut[i],150,"M_{TOF}^{2}_%s_2#sigma_TPCcut_DCAxyCut;M_{TOF}^{2}/Z^{2} (GeV^{2}/c^{4});p/|Z| (GeV/c)",name[i]);
214   }
215
216   for (Int_t i=0;i<18;i++) fM2vsP_DCAxyCut[i] = new TH2F(name_par_MvsP_DCAxyCut[i],name_title_MvsP_DCAxyCut[i],8000,0,10,200,0,10);
217
218
219   Char_t par_name_nbin[nbin][30];
220   
221   snprintf(par_name_nbin[0],30,"0.4<Pt<0.5");
222   snprintf(par_name_nbin[1],30,"0.5<Pt<0.6");
223   snprintf(par_name_nbin[2],30,"0.6<Pt<0.7");
224   snprintf(par_name_nbin[3],30,"0.7<Pt<0.8");
225   snprintf(par_name_nbin[4],30,"0.8<Pt<0.9");
226   snprintf(par_name_nbin[5],30,"0.9<Pt<1.0");
227   snprintf(par_name_nbin[6],30,"1.0<Pt<1.1");
228   snprintf(par_name_nbin[7],30,"1.1<Pt<1.2");
229   snprintf(par_name_nbin[8],30,"1.2<Pt<1.3");
230   snprintf(par_name_nbin[9],30,"1.3<Pt<1.4");
231   snprintf(par_name_nbin[10],30,"1.4<Pt<1.5");
232   snprintf(par_name_nbin[11],30,"1.5<Pt<1.6");
233   snprintf(par_name_nbin[12],30,"1.6<Pt<1.7");
234   snprintf(par_name_nbin[13],30,"1.7<Pt<1.8");
235   snprintf(par_name_nbin[14],30,"1.8<Pt<1.9");
236   snprintf(par_name_nbin[15],30,"1.9<Pt<2.0");
237   snprintf(par_name_nbin[16],30,"2.0<Pt<2.1");
238   snprintf(par_name_nbin[17],30,"2.1<Pt<2.2");
239   snprintf(par_name_nbin[18],30,"2.2<Pt<2.3");
240   snprintf(par_name_nbin[19],30,"2.3<Pt<2.4");
241   snprintf(par_name_nbin[20],30,"2.4<Pt<2.5");
242   snprintf(par_name_nbin[21],30,"2.5<Pt<2.6");
243   snprintf(par_name_nbin[22],30,"2.6<Pt<2.7");
244   snprintf(par_name_nbin[23],30,"2.7<Pt<2.8");
245   snprintf(par_name_nbin[24],30,"2.8<Pt<2.9");
246   snprintf(par_name_nbin[25],30,"2.9<Pt<3.0");
247   snprintf(par_name_nbin[26],30,"3.0<Pt<3.1");
248   snprintf(par_name_nbin[27],30,"3.1<Pt<3.2");
249   snprintf(par_name_nbin[28],30,"3.2<Pt<3.3");
250   snprintf(par_name_nbin[29],30,"3.3<Pt<3.4");
251   snprintf(par_name_nbin[30],30,"3.4<Pt<3.5");
252   snprintf(par_name_nbin[31],30,"3.5<Pt<3.6");
253   snprintf(par_name_nbin[32],30,"3.6<Pt<3.7");
254   snprintf(par_name_nbin[33],30,"3.7<Pt<3.8");
255   snprintf(par_name_nbin[34],30,"3.8<Pt<3.9");
256   snprintf(par_name_nbin[35],30,"3.9<Pt<4.0");
257   snprintf(par_name_nbin[36],30,"4.0<Pt<4.1");
258   snprintf(par_name_nbin[37],30,"4.1<Pt<4.2");
259   snprintf(par_name_nbin[38],30,"4.2<Pt<4.3");
260   snprintf(par_name_nbin[39],30,"4.3<Pt<4.4");
261   snprintf(par_name_nbin[40],30,"4.4<Pt<4.5");
262   snprintf(par_name_nbin[41],30,"4.5<Pt<4.6");
263   snprintf(par_name_nbin[42],30,"4.6<Pt<4.7");
264   snprintf(par_name_nbin[43],30,"4.7<Pt<4.8");
265   snprintf(par_name_nbin[44],30,"4.8<Pt<4.9");
266   snprintf(par_name_nbin[45],30,"4.9<Pt<5.0");
267
268   
269   Char_t nameDCAxy[18][nbin][120];
270   Char_t titleDCAxy[18][nbin][120];
271
272   Char_t nameM2CutDCAxy[18][nbin][120];
273   Char_t titleM2CutDCAxy[18][nbin][120];
274  
275   Char_t nameM2CutGroundDCAxy[18][nbin][120];
276   Char_t titleM2CutGroundDCAxy[18][nbin][120];
277
278   
279   for(Int_t iS=0;iS<18;iS++) {
280     for(Int_t j=0;j<nbin;j++) {
281       snprintf(nameDCAxy[iS][j],120,"hDCAxy_%s_%s",name[iS],par_name_nbin[j]);
282       snprintf(titleDCAxy[iS][j],120,"hDCAxy_%s_%s;DCA_{xy} (cm)",name[iS],par_name_nbin[j]);
283     
284       snprintf(nameM2CutDCAxy[iS][j],120,"hM2_CutDCAxy_%s_%s",name[iS],par_name_nbin[j]);
285       snprintf(titleM2CutDCAxy[iS][j],120,"hM2_CutDCAxy_%s_%s;M^{2}_{TOF} (GeV^{2}/c^{4})",name[iS],par_name_nbin[j]);
286       
287       snprintf(nameM2CutGroundDCAxy[iS][j],120,"hM2_GroundCatDCAxy_%s_%s",name[iS],par_name_nbin[j]);
288       snprintf(titleM2CutGroundDCAxy[iS][j],120,"hM2_GroundCatDCAxy_%s_%s;M^{2}_{TOF} (GeV^{2}/c^{4})",name[iS],par_name_nbin[j]);
289     }
290   }
291  
292   for(Int_t iS=0;iS<18;iS++) {
293     for(Int_t j=0;j<nbin;j++) {
294       hDCAxy[iS][j] = new TH1D(nameDCAxy[iS][j],titleDCAxy[iS][j],1000,-2.5,2.5);//125 bins
295       hDCAxy[iS][j]->GetXaxis()->CenterTitle();
296       //hM2CutDCAxy[iS][j] = new TH1D(nameM2CutDCAxy[iS][j],titleM2CutDCAxy[iS][j],4000,0,10);
297       //hM2CutGroundDCAxy[iS][j] = new TH1D(nameM2CutGroundDCAxy[iS][j],titleM2CutGroundDCAxy[iS][j],4000,0,10);
298     }
299   }
300   
301   //for e,#mu,#pi and antiparticle (e and #mu will not be drawn)
302   //the binning is chosen for #pi distribution:
303   for(Int_t iSp=0;iSp<3;iSp++) {
304     for(Int_t j=0;j<nbin;j++) {
305       hM2CutDCAxy[iSp][j] = new TH1D(nameM2CutDCAxy[iSp][j],titleM2CutDCAxy[iSp][j],600,-0.1,0.5);
306       hM2CutGroundDCAxy[iSp][j] = new TH1D(nameM2CutGroundDCAxy[iSp][j],titleM2CutGroundDCAxy[iSp][j],600,-0.1,0.5);
307       hM2CutDCAxy[iSp+9][j] = new TH1D(nameM2CutDCAxy[iSp+9][j],titleM2CutDCAxy[iSp+9][j],600,-0.1,0.5);
308       hM2CutGroundDCAxy[iSp+9][j] = new TH1D(nameM2CutGroundDCAxy[iSp+9][j],titleM2CutGroundDCAxy[iSp+9][j],600,-0.1,0.5);
309     }
310   }
311
312   for(Int_t j=0;j<nbin;j++) {
313     hM2CutDCAxy[3][j] = new TH1D(nameM2CutDCAxy[3][j],titleM2CutDCAxy[3][j],400,0,1);
314     hM2CutGroundDCAxy[3][j] = new TH1D(nameM2CutGroundDCAxy[3][j],titleM2CutGroundDCAxy[3][j],400,0,1);
315     hM2CutDCAxy[3+9][j] = new TH1D(nameM2CutDCAxy[3+9][j],titleM2CutDCAxy[3+9][j],400,0,1);
316     hM2CutGroundDCAxy[3+9][j] = new TH1D(nameM2CutGroundDCAxy[3+9][j],titleM2CutGroundDCAxy[3+9][j],400,0,1);
317   }
318
319   for(Int_t j=0;j<nbin;j++) {
320     hM2CutDCAxy[4][j] = new TH1D(nameM2CutDCAxy[4][j],titleM2CutDCAxy[4][j],500,0,4);
321     hM2CutGroundDCAxy[4][j] = new TH1D(nameM2CutGroundDCAxy[4][j],titleM2CutGroundDCAxy[4][j],500,0,4);
322     hM2CutDCAxy[4+9][j] = new TH1D(nameM2CutDCAxy[4+9][j],titleM2CutDCAxy[4+9][j],500,0,4);
323     hM2CutGroundDCAxy[4+9][j] = new TH1D(nameM2CutGroundDCAxy[4+9][j],titleM2CutGroundDCAxy[4+9][j],500,0,4);
324   }
325
326   for(Int_t j=0;j<nbin;j++) {
327     hM2CutDCAxy[5][j] = new TH1D(nameM2CutDCAxy[5][j],titleM2CutDCAxy[5][j],500,0,6);
328     hM2CutGroundDCAxy[5][j] = new TH1D(nameM2CutGroundDCAxy[5][j],titleM2CutGroundDCAxy[5][j],500,0,6);
329     hM2CutDCAxy[5+9][j] = new TH1D(nameM2CutDCAxy[5+9][j],titleM2CutDCAxy[5+9][j],500,0,6);
330     hM2CutGroundDCAxy[5+9][j] = new TH1D(nameM2CutGroundDCAxy[5+9][j],titleM2CutGroundDCAxy[5+9][j],500,0,6);
331   }
332
333   for(Int_t j=0;j<nbin;j++) {
334     hM2CutDCAxy[6][j] = new TH1D(nameM2CutDCAxy[6][j],titleM2CutDCAxy[6][j],1000,0,12);
335     hM2CutGroundDCAxy[6][j] = new TH1D(nameM2CutGroundDCAxy[6][j],titleM2CutGroundDCAxy[6][j],1000,0,12);
336     hM2CutDCAxy[6+9][j] = new TH1D(nameM2CutDCAxy[6+9][j],titleM2CutDCAxy[6+9][j],1000,0,12);
337     hM2CutGroundDCAxy[6+9][j] = new TH1D(nameM2CutGroundDCAxy[6+9][j],titleM2CutGroundDCAxy[6+9][j],1000,0,12);
338   }
339
340   for(Int_t j=0;j<nbin;j++) {
341     hM2CutDCAxy[7][j] = new TH1D(nameM2CutDCAxy[7][j],titleM2CutDCAxy[7][j],200,0,4);
342     hM2CutGroundDCAxy[7][j] = new TH1D(nameM2CutGroundDCAxy[7][j],titleM2CutGroundDCAxy[7][j],200,0,4);
343     hM2CutDCAxy[7+9][j] = new TH1D(nameM2CutDCAxy[7+9][j],titleM2CutDCAxy[7+9][j],200,0,4);
344     hM2CutGroundDCAxy[7+9][j] = new TH1D(nameM2CutGroundDCAxy[7+9][j],titleM2CutGroundDCAxy[7+9][j],200,0,4);
345   }
346
347   for(Int_t j=0;j<nbin;j++) {
348     hM2CutDCAxy[8][j] = new TH1D(nameM2CutDCAxy[8][j],titleM2CutDCAxy[8][j],600,0,6);
349     hM2CutGroundDCAxy[8][j] = new TH1D(nameM2CutGroundDCAxy[8][j],titleM2CutGroundDCAxy[8][j],600,0,6);
350     hM2CutDCAxy[8+9][j] = new TH1D(nameM2CutDCAxy[8+9][j],titleM2CutDCAxy[8+9][j],600,0,6);
351     hM2CutGroundDCAxy[8+9][j] = new TH1D(nameM2CutGroundDCAxy[8+9][j],titleM2CutGroundDCAxy[8+9][j],600,0,6);
352   }
353
354  
355   fList1->Add(hNevent);
356   fList1->Add(hZvertex);
357   fList1->Add(hTOFSignalPion);
358   for(Int_t iS=0;iS<18;iS++) fList1->Add(fNsigmaTPCvsP_kTOFtrue[iS]);
359    
360   for(Int_t i=0;i<3;i++) fList1->Add(fdEdxVSp[i]);
361   for(Int_t i=0;i<9;i++) fList1->Add(hDeDxExp[i]);
362   fList1->Add(fBetaTofVSp);
363   for(Int_t i=0;i<9;i++) fList1->Add(hBetaExp[i]);
364   for(Int_t i=0;i<9;i++) fList1->Add(fNsigmaTPC[i]);
365   for(Int_t i=0;i<9;i++) fList1->Add(fNsigmaTOF[i]);
366
367   for(Int_t i=0;i<10;i++) fList1->Add(fM2vsZ[i]);
368   
369   for(Int_t i=0;i<3;i++) fList1->Add(fM2vsP_NoTpcCut[i]);
370   for(Int_t i=0;i<18;i++) fList1->Add(fM2vsP[i]);
371   
372   for(Int_t i=0;i<3;i++) fList1->Add(fM2vsP_NoTpcCut_DCAxyCut[i]);
373   for(Int_t i=0;i<18;i++) fList1->Add(fM2vsP_DCAxyCut[i]);
374
375   /*
376   for(Int_t j=0;j<nbin;j++) {//electron
377     fList1->Add(hDCAxy[0][j]);
378     fList1->Add(hM2CutDCAxy[0][j]);
379     fList1->Add(hM2CutGroundDCAxy[0][j]);
380     fList1->Add(hDCAxy[9][j]);
381     fList1->Add(hM2CutDCAxy[9][j]);
382     fList1->Add(hM2CutGroundDCAxy[9][j]);
383   }
384
385   for(Int_t j=0;j<nbin;j++) {//muon
386     fList1->Add(hDCAxy[1][j]);
387     fList1->Add(hM2CutDCAxy[1][j]);
388     fList1->Add(hM2CutGroundDCAxy[1][j]);
389     fList1->Add(hDCAxy[10][j]);
390     fList1->Add(hM2CutDCAxy[10][j]);
391     fList1->Add(hM2CutGroundDCAxy[10][j]);
392   }
393   */
394
395   for(Int_t j=0;j<nbin;j++) {//pion
396     fList1->Add(hDCAxy[2][j]);
397     fList1->Add(hM2CutDCAxy[2][j]);
398     fList1->Add(hM2CutGroundDCAxy[2][j]);
399     fList1->Add(hDCAxy[11][j]);
400     fList1->Add(hM2CutDCAxy[11][j]);
401     fList1->Add(hM2CutGroundDCAxy[11][j]);
402   }
403
404   for(Int_t j=0;j<nbin;j++) {//kaon
405     fList1->Add(hDCAxy[3][j]);
406     fList1->Add(hM2CutDCAxy[3][j]);
407     fList1->Add(hM2CutGroundDCAxy[3][j]);
408     fList1->Add(hDCAxy[12][j]);
409     fList1->Add(hM2CutDCAxy[12][j]);
410     fList1->Add(hM2CutGroundDCAxy[12][j]);
411   }
412
413   for(Int_t j=0;j<nbin;j++) {//proton
414     fList1->Add(hDCAxy[4][j]);
415     fList1->Add(hM2CutDCAxy[4][j]);
416     fList1->Add(hM2CutGroundDCAxy[4][j]);
417     fList1->Add(hDCAxy[13][j]);
418     fList1->Add(hM2CutDCAxy[13][j]);
419     fList1->Add(hM2CutGroundDCAxy[13][j]);
420   }
421
422   for(Int_t j=0;j<nbin;j++) {//deuteron
423     fList1->Add(hDCAxy[5][j]);
424     fList1->Add(hM2CutDCAxy[5][j]);
425     fList1->Add(hM2CutGroundDCAxy[5][j]);
426     fList1->Add(hDCAxy[14][j]);
427     fList1->Add(hM2CutDCAxy[14][j]);
428     fList1->Add(hM2CutGroundDCAxy[14][j]);
429   }
430
431   /*
432   for(Int_t j=0;j<nbin;j++) {//triton
433     fList1->Add(hDCAxy[6][j]);
434     fList1->Add(hM2CutDCAxy[6][j]);
435     fList1->Add(hM2CutGroundDCAxy[6][j]);
436     fList1->Add(hDCAxy[15][j]);
437     fList1->Add(hM2CutDCAxy[15][j]);
438     fList1->Add(hM2CutGroundDCAxy[15][j]);
439   }
440   */
441
442   for(Int_t j=0;j<nbin;j++) {//He3
443     fList1->Add(hDCAxy[7][j]);
444     fList1->Add(hM2CutDCAxy[7][j]);
445     fList1->Add(hM2CutGroundDCAxy[7][j]);
446     fList1->Add(hDCAxy[16][j]);
447     fList1->Add(hM2CutDCAxy[16][j]);
448     fList1->Add(hM2CutGroundDCAxy[16][j]);
449   }
450   
451   /*
452   for(Int_t j=0;j<nbin;j++) {//4He
453     fList1->Add(hDCAxy[8][j]);
454     fList1->Add(hM2CutDCAxy[8][j]);
455     fList1->Add(hM2CutGroundDCAxy[8][j]);
456     fList1->Add(hDCAxy[17][j]);
457     fList1->Add(hM2CutDCAxy[17][j]);
458     fList1->Add(hM2CutGroundDCAxy[17][j]);
459   }
460   */
461
462   //fTrackFilter = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
463   
464   // Post output data.
465   PostData(1, fList1);
466 }
467 //______________________________________________________________________________
468 void AliAnalysisNucleiMass::UserExec(Option_t *) 
469 {
470   // Main loop
471   // Called for each event
472   
473   fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
474   fESD = dynamic_cast<AliESDEvent*>(InputEvent());
475   if(!fAOD && !fESD){
476     Printf("%s:%d AODEvent and ESDEvent not found in Input Manager",(char*)__FILE__,__LINE__);
477     this->Dump();
478     return;
479   }
480   
481   if(fESD) fEvent = fESD;
482   else fEvent = fAOD;
483
484   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
485   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
486   fPIDResponse=inputHandler->GetPIDResponse(); // data member di tipo "const AliPIDResponse *fPIDResponse;"
487
488   //Centrality
489   Float_t v0Centr  = -10.;
490   Float_t trkCentr  = -10.;
491   AliCentrality *centrality = fEvent->GetCentrality();
492   if (centrality){
493     v0Centr  = centrality->GetCentralityPercentile("V0M"); // VZERO
494     trkCentr = centrality->GetCentralityPercentile("TRK"); // TPC
495   }
496
497   const AliAODVertex* vtxAOD = fAOD->GetPrimaryVertex();
498   
499   Float_t zvtx = 10000.0;
500
501   if(vtxAOD->GetNContributors()>0)
502     zvtx = vtxAOD->GetZ();
503   
504   if(TMath::Abs(v0Centr - trkCentr) < 5.0 && TMath::Abs(zvtx) < 10.0){ // consistency cut on centrality selection AND d(zPrimaryVertez;NominalPointInteraction)<10cm
505     
506     //Bool_t isTrack=1;
507     
508     Int_t nTracks = fEvent->GetNumberOfTracks();
509     
510     if(v0Centr>=fCentrality[0] && v0Centr<=fCentrality[1]) {//window cut centrality open
511       
512       hZvertex->Fill(zvtx);
513       hNevent->Fill(v0Centr);
514       
515       for(Int_t iT = 0; iT < nTracks; iT++) { // loop on the tracks
516         AliVTrack* track = (AliVTrack *) fEvent->GetTrack(iT);
517         
518         if (!track){
519           track->Delete();
520           //isTrack=0;
521           continue;
522         }
523         
524         Bool_t trkFlag = 0;
525         
526         if(fAOD)
527           trkFlag = ((AliAODTrack *) track)->TestFilterBit(FilterBit);
528         
529           //TestFilterBit(16) -- Standard Cuts with very loose DCA: GetStandardITSTPCTrackCuts2011(kFALSE) && SetMaxDCAToVertexXY(2.4) && SetMaxDCAToVertexZ(3.2) && SetDCaToVertex2D(kTRUE)
530         
531           //TestFilterBit(32) (STARDARD) -- Standard Cuts with very tight DCA cut ( 7sigma^primaries: 7*(0.0015+0.0050/pt^1.1) ) : GetStandardITSTPCTrackCuts2011(). 
532
533         else{
534           //trkFlag = fTrackFilter->IsSelected(((AliESDtrack *) track));
535         }
536         
537         Int_t NTpcCls=track->GetTPCNcls();
538         if(NTpcCls>NminTPCcluster) kTPC=kTRUE;
539         else kTPC=kFALSE;
540
541         if ((TMath::Abs(track->Eta()) > 0.8) || (track->Pt() < 0.2) || !trkFlag || !kTPC){
542           continue;
543         }
544         
545         Double_t b[2] = {-99., -99.};
546         Double_t bCov[3] = {-99., -99., -99.};
547         if (!track->PropagateToDCA(fEvent->GetPrimaryVertex(), fEvent->GetMagneticField(), 100., b, bCov))
548           continue;
549         
550 //      Float_t Eta = TMath::Abs(track->Eta());
551         Float_t charge = (Float_t)track->Charge();
552         Float_t p = track->P();
553         Float_t pt = track->Pt();
554         Float_t dedx = track->GetTPCsignal();
555         Float_t tof  = track->GetTOFsignal()-fPIDResponse->GetTOFResponse().GetStartTime(p);
556         Float_t pTPC = track->GetTPCmomentum();
557         Float_t beta = 0.0;
558         Float_t M2 = 1000.0;
559         Float_t M = 1000.0;
560         Float_t Z2 = 1000.0;
561         Float_t DCAxy = b[1];
562         Float_t DCAz = b[0];
563
564         if(TMath::Abs(DCAz)>DCAzCUT)//CUT ON DCAz
565           continue;
566         
567         Bool_t kTpcPure;
568         kTpcPure = track->GetTPCsignal()>10;
569         
570         kTOF = (track->GetStatus() & AliVTrack::kTOFout) && (track->GetStatus() & AliVTrack::kTIME);
571
572         Float_t nsigmaTPC[9];
573         Float_t nsigmaTOF[9];
574         
575         Double_t expdedx[9];
576 //      Double_t dedxRes[9];
577         
578         for(Int_t iS=0;iS < 9;iS++){ //TPC expected signal
579           //expdedx[iS] = fPIDResponse->GetTPCResponse().GetExpectedSignal(pTPC, (AliPID::EParticleType) iS);// Deprecated function: Temporary solution to measure also He (The other method has a small bag) . This must be replace from:
580           expdedx[iS] = fPIDResponse->GetTPCResponse().GetExpectedSignal(track, (AliPID::EParticleType) iS, AliTPCPIDResponse::kdEdxDefault, kTRUE);
581 //        dedxRes[iS] = fPIDResponse->GetTPCResponse().GetExpectedSigma(track, (AliPID::EParticleType) iS, AliTPCPIDResponse::kdEdxDefault, kTRUE);
582           //if(iS>6) dedxRes[iS] = expdedx[iS]*0.08; //Temporary solution to measure also He. This line must be removed 
583         }
584         
585         if(kTpcPure) {
586           for(Int_t iS=0;iS < 9;iS++){
587             nsigmaTPC[iS] = fPIDResponse->NumberOfSigmasTPC(track,(AliPID::EParticleType) iS);
588             //if(iS>6) nsigmaTPC[iS] = (dedx-expdedx[iS])/dedxRes[iS]; // Temporary solution to measure also He. This line must be removed.
589             fNsigmaTPC[iS]->Fill(pt,nsigmaTPC[iS]);
590             hDeDxExp[iS]->Fill(pTPC,expdedx[iS]);
591           }
592           fdEdxVSp[0]->Fill(pTPC,dedx);
593           if(charge>0) fdEdxVSp[1]->Fill(pTPC,dedx);
594           else fdEdxVSp[2]->Fill(pTPC,dedx);
595         }
596         
597         Float_t massOverZ[9] = {0.000511,0.105658,0.139570,0.493677,0.938272,1.877837,2.817402,1.408701,1.877837};
598         
599         Double_t exptimes[9]; // TOF expected times
600         track->GetIntegratedTimes(exptimes);
601         exptimes[5] = exptimes[0] / p * massOverZ[5] * TMath::Sqrt(1+p*p/massOverZ[5]/massOverZ[5]);
602         exptimes[6] = exptimes[0] / p * massOverZ[6] * TMath::Sqrt(1+p*p/massOverZ[6]/massOverZ[6]);
603         exptimes[7] = exptimes[0] / p * massOverZ[7] * TMath::Sqrt(1+p*p/massOverZ[7]/massOverZ[7]);
604         exptimes[8] = exptimes[0] / p * massOverZ[8] * TMath::Sqrt(1+p*p/massOverZ[8]/massOverZ[8]);
605         
606 //      Double_t tofRes[9] = {1000,1000,1000,1000,1000,1000,1000,1000,1000};
607 //      if(kTOF){
608 //        for(Int_t iS=0;iS < 9;iS++) tofRes[iS] = fPIDResponse->GetTOFResponse().GetExpectedSigma(p, exptimes[iS], massOverZ[iS]);
609 //      }
610         
611         beta=exptimes[0];
612         
613         if(kTOF) {
614           for(Int_t iS=0;iS < 9;iS++){
615             nsigmaTOF[iS] = fPIDResponse->NumberOfSigmasTOF(track,(AliPID::EParticleType) iS);
616             fNsigmaTOF[iS]->Fill(pt,nsigmaTOF[iS]);
617             hBetaExp[iS]->Fill(p,beta/exptimes[iS]);
618           }
619           if(pt>0.9 && pt<1.0) hTOFSignalPion->Fill(tof-exptimes[2]);
620           beta=beta/tof;
621           fBetaTofVSp->Fill(p,beta);
622         }
623         
624         Int_t stdFlagPid[9] = {1,2,4,8,16,32,64,128,256};//e,#mu,#pi,K,p,d,t,3He,4He
625         Int_t FlagPid = 0;
626         Float_t binPt[nbin+1];
627
628         for(Int_t i=0;i<nbin+1;i++) {
629           binPt[i]=0.4+i*0.1;
630         }
631         
632         Float_t binCutPt[10] = {0.3,0.5,1.0,1.5,2.0,2.5,3.0,3.5,4.0,4.5};
633
634         if(kTOF && kTpcPure) {
635           
636           M2 = (p*p*(1-beta*beta))/(beta*beta);
637           
638           fM2vsP_NoTpcCut[0]->Fill(M2,p);
639           if(TMath::Abs(DCAxy)<DCAxyCUT) fM2vsP_NoTpcCut_DCAxyCut[0]->Fill(M2,p);
640           
641           if(M2>0.0) {
642             M=TMath::Sqrt(M2);
643             Z2 = TMath::Power(dedx/fPIDResponse->GetTPCResponse().GetExpectedSignal(pTPC*massOverZ[4]/M, AliPID::kProton),0.862);
644             fM2vsZ[0]->Fill(charge*TMath::Sqrt(Z2),M2);
645             
646             for(Int_t i=0;i<9;i++) {
647               if(pt>binCutPt[i] && pt<binCutPt[i+1]){
648                 fM2vsZ[i+1]->Fill(charge*TMath::Sqrt(Z2),M2);
649                 break;
650               }
651             }
652           }
653           
654           if(charge>0) {
655             fM2vsP_NoTpcCut[1]->Fill(M2,p);
656             if(TMath::Abs(DCAxy)<DCAxyCUT) fM2vsP_NoTpcCut_DCAxyCut[1]->Fill(M2,p);
657             for(Int_t iS=0;iS < 9;iS++){
658               fNsigmaTPCvsP_kTOFtrue[iS]->Fill(p,nsigmaTPC[iS]);
659             }
660           }
661           else {//else charge<0
662             fM2vsP_NoTpcCut[2]->Fill(M2,p);
663             if(TMath::Abs(DCAxy)<DCAxyCUT) fM2vsP_NoTpcCut_DCAxyCut[2]->Fill(M2,p);
664             
665             for(Int_t iS=0;iS < 9;iS++){
666               fNsigmaTPCvsP_kTOFtrue[iS+9]->Fill(p,nsigmaTPC[iS]);
667             }
668           }
669           
670           for(Int_t iS=0;iS<9;iS++) {
671             if(TMath::Abs(nsigmaTPC[iS])<2) {
672               FlagPid += ((Int_t)TMath::Power(2,iS));
673             }
674           }
675             
676           for(Int_t iS=0;iS<9;iS++) {
677             if(FlagPid & stdFlagPid[iS] || !kTPCcut) {
678               if(charge>0) {
679                 fM2vsP[iS]->Fill(M2,p);
680                 for(Int_t j=0;j<nbin;j++) {
681                   if(pt>binPt[j] && pt<binPt[j+1]) {
682                     hDCAxy[iS][j]->Fill(DCAxy);
683                     hDCAxy[iS][j]->Fill(-DCAxy);
684                     if(TMath::Abs(DCAxy)<DCAxyCUT) {
685                       hM2CutDCAxy[iS][j]->Fill(M2);
686                       fM2vsP_DCAxyCut[iS]->Fill(M2,p);
687                     }
688                     if(TMath::Abs(DCAxy+0.5)<DCAxyCUT) hM2CutGroundDCAxy[iS][j]->Fill(M2);
689                     break;
690                   }
691                 }
692               }
693               else {//if(charge<0)
694                 fM2vsP[iS+9]->Fill(M2,p);
695                 for(Int_t j=0;j<nbin;j++) {
696                   if(pt>binPt[j] && pt<binPt[j+1]) {
697                     hDCAxy[iS+9][j]->Fill(DCAxy);
698                     hDCAxy[iS+9][j]->Fill(-DCAxy);
699                     if(TMath::Abs(DCAxy)<DCAxyCUT) {
700                       hM2CutDCAxy[iS+9][j]->Fill(M2);
701                       fM2vsP_DCAxyCut[iS]->Fill(M2,p);
702                     }
703                     if(TMath::Abs(DCAxy+0.5)<DCAxyCUT) hM2CutGroundDCAxy[iS+9][j]->Fill(M2);
704                     break;
705                   }
706                 }
707               }
708             }
709           }
710         }// close (KTOF && kTpcPure request)
711         
712         
713       } // end track loop
714       
715     } //window cut centrality close
716     
717   }  
718   
719 }
720
721 //_____________________________________________________________________________
722 void AliAnalysisNucleiMass::Terminate(Option_t *)
723
724   // Terminate loop
725   Printf("Terminate()");
726 }