]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/macros/AddTaskCFVertexingHF3ProngLc.C
Centrality bins for 2011 PbPb run (Davide)
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / macros / AddTaskCFVertexingHF3ProngLc.C
1 //DEFINITION OF A FEW CONSTANTS
2 const Double_t ymin  = -1.2 ;
3 const Double_t ymax  =  1.2 ;
4 const Double_t cosmin = -0.7;
5 const Double_t cosmax =  1.05;
6 const Double_t cTmin = 0;  // micron
7 const Double_t cTmax = 500;  // micron
8 const Double_t phimin = 0.0;  
9 const Int_t    mintrackrefsTPC = 2 ;
10 const Int_t    mintrackrefsITS = 3 ;
11 const Int_t    charge  = 1 ;
12 const Int_t    minclustersTPC = 50 ;
13 // cuts
14 const Double_t ptmin = 0.1;
15 const Double_t ptmax = 9999.;
16 const Double_t etamin = -0.9;
17 const Double_t etamax = 0.9;
18 const Double_t zvtxmin = -15;
19 const Double_t zvtxmax = 15;
20 const Int_t    minITSClusters = 5;
21
22 const Float_t centmin_0_10 = 0.;
23 const Float_t centmax_0_10 = 10.;
24 const Float_t centmin_10_60 = 10.;
25 const Float_t centmax_10_60 = 60.;
26 const Float_t centmin_60_100 = 60.;
27 const Float_t centmax_60_100 = 100.;
28 const Float_t centmax = 100.;
29 const Float_t fakemin = -0.5;
30 const Float_t fakemax = 2.5.;
31 const Double_t distTwoPartmin=0;
32 const Double_t distTwoPartmax=600;
33 const Double_t dispVtxmin = 0;
34 const Double_t dispVtxmax = 600;
35 const Double_t sumd02min = 0.;
36 const Double_t sumd02max = 50000.;
37 const Float_t cosminXY = 0.95;
38 const Float_t cosmaxXY = 1.0;
39 const Float_t normDecLXYmin = 0;
40 const Float_t normDecLXYmax = 20;
41 const Float_t multmin_0_20 = 0;
42 const Float_t multmax_0_20 = 20;
43 const Float_t multmin_20_50 = 20;
44 const Float_t multmax_20_50 = 50;
45 const Float_t multmin_50_102 = 50;
46 const Float_t multmax_50_102 = 102;
47
48
49 //----------------------------------------------------
50
51 //AliCFTaskVertexingHF *AddTaskCFVertexingHF3ProngLc(const char* cutFile = "./cuts4LctopKpi.root", Int_t configuration = AliCFTaskVertexingHF::kSnail, Bool_t isKeepDfromB=kFALSE, Bool_t isKeepDfromBOnly=kFALSE, Int_t pdgCode = 4122, Char_t isSign = 2)
52 AliCFContainer *AddTaskCFVertexingHF3ProngLc(const char* cutFile = "./cuts4LctopKpi.root", Int_t configuration = AliCFTaskVertexingHF::kSnail, Bool_t isKeepDfromB=kFALSE, Bool_t isKeepDfromBOnly=kFALSE, Int_t pdgCode = 4122, Char_t isSign = 2)
53 {
54         printf("Addig CF task using cuts from file %s\n",cutFile);
55         if (configuration == AliCFTaskVertexingHF::kSnail){
56                 printf("The configuration is set to be SLOW --> all the variables will be used to fill the CF\n");
57         }
58         else if (configuration == AliCFTaskVertexingHF::kCheetah){
59                 printf("The configuration is set to be FAST --> using only pt, y, ct, phi, zvtx, centrality, fake, multiplicity to fill the CF\n");
60         }
61         else{
62                 printf("The configuration is not defined! returning\n");
63                 return;
64         }
65         
66         gSystem->Sleep(2000);
67         
68         // isSign = 0 --> D0 only
69         // isSign = 1 --> D0bar only
70         // isSign = 2 --> D0 + D0bar
71         
72         TString expected;
73         if (isSign == 0 && pdgCode < 0){
74                 AliError(Form("Error setting PDG code (%d) and sign (0 --> particle (%d) only): they are not compatible, returning",pdgCode));
75                 return 0x0;
76         }
77         else if (isSign == 1 && pdgCode > 0){
78                 AliError(Form("Error setting PDG code (%d) and sign (1 --> antiparticle (%d) only): they are not compatible, returning",pdgCode));
79                 return 0x0;
80         }
81         else if (isSign > 2 || isSign < 0){
82                 AliError(Form("Sign not valid (%d, possible values are 0, 1, 2), returning"));
83                 return 0x0;
84         }
85         
86         TFile* fileCuts = TFile::Open(cutFile);
87         if(!fileCuts || (fileCuts && !fileCuts->IsOpen())){ 
88           AliError("Wrong cut file");
89           return 0x0;
90         }
91
92         AliRDHFCutsLctopKpi *cutsLctopKpi = (AliRDHFCutsLctopKpi*)fileCuts->Get("LctopKpiProdCuts");
93         
94         // check that the fKeepD0fromB flag is set to true when the fKeepD0fromBOnly flag is true
95         //  for now the binning is the same than for all D's
96         if(isKeepDfromBOnly) isKeepDfromB = true;
97         
98         /*
99           Double_t ptmin_0_4;
100           Double_t ptmax_0_4;
101           Double_t ptmin_4_8;
102           Double_t ptmax_4_8;
103           Double_t ptmin_8_10;
104           Double_t ptmax_8_10;
105           
106           if(!isKeepDfromB){
107           ptmin_0_4 =  0.0 ;
108           ptmax_0_4 =  4.0 ;
109           ptmin_4_8 =  4.0 ;
110           ptmax_4_8 =  8.0 ;
111           ptmin_8_10 =  8.0 ;
112           ptmax_8_10 =  10.0 ;
113           } else{
114           ptmin_0_4 =  0.0 ;
115           ptmax_0_4 =  3.0 ;
116           ptmin_4_8 =  3.0 ;
117           ptmax_4_8 =  5.0 ;
118           ptmin_8_10 =  5.0 ;
119           ptmax_8_10 =  10.0 ;
120           }
121         */
122
123         //CONTAINER DEFINITION
124         Info("AliCFTaskVertexingHF","SETUP CONTAINER");
125
126         const Double_t phimax = 2*TMath::Pi();
127
128         //Setting up the container grid... 
129         UInt_t nstep = 10; //number of selection steps: MC with limited acceptance, MC, Acceptance, Vertex, Refit, Reco (no cuts), RecoAcceptance, RecoITSClusters (RecoAcceptance included), RecoPPR (RecoAcceptance+RecoITSCluster included), RecoPID 
130 //      const Int_t nbinpt_0_4  = 8 ; //bins in pt from 0 to 4 GeV
131 //      const Int_t nbinpt_4_8  = 4 ; //bins in pt from 4 to 8 GeV
132 //      const Int_t nbinpt_8_10  = 1 ; //bins in pt from 8 to 10 GeV
133
134 /*
135         Int_t nbinpt_0_4;
136         Int_t nbinpt_4_8;
137         Int_t nbinpt_8_10;
138         if (!isKeepDfromB){
139           nbinpt_0_4  = 8 ; //bins in pt from 0 to 4 GeV
140           nbinpt_4_8  = 4 ; //bins in pt from 4 to 8 GeV
141           nbinpt_8_10  = 1 ; //bins in pt from 8 to 10 GeV
142         }else{
143           nbinpt_0_4  = 3 ; //bins in pt from 0 to 3 GeV
144           nbinpt_4_8  = 1 ; //bins in pt from 3 to 5 GeV
145           nbinpt_8_10  = 1 ; //bins in pt from 5 to 10 GeV
146         }
147 */
148         const Int_t nbinpt = cutsLctopKpi->GetNPtBins(); // bins in pT
149         printf("pT: nbin (from cuts file) = %d\n",nbinpt);
150         const Int_t nbiny  = 24 ; //bins in y
151         const Int_t nbinphi  = 18 ; //bins in phi
152         const Int_t nbincT  = 25 ; //bins in cT 
153         const Int_t nbinpointing  = 350 ; //bins in cosPointingAngle    
154         const Int_t nbinpTpi_0_4  = 8 ; //bins in ptPi from 0 to 4 GeV
155         const Int_t nbinpTpi_4_8  = 4 ; //bins in ptPi from 4 to 8 GeV
156         const Int_t nbinpTpi_8_10  = 1 ; //bins in ptPi from 8 to 10 GeV
157         const Int_t nbinpTk_0_4  = 8 ; //bins in ptKa from 0 to 4 GeV
158         const Int_t nbinpTk_4_8  = 4 ; //bins in ptKa from 4 to 8 GeV
159         const Int_t nbinpTk_8_10  = 1 ; //bins in ptKa from 8 to 10 GeV
160         const Int_t nbinpTpi2_0_4  = 8 ; //bins in ptpi2 from 0 to 4 GeV
161         const Int_t nbinpTpi2_4_8  = 4 ; //bins in ptpi2 from 4 to 8 GeV
162         const Int_t nbinpTpi2_8_10  = 1 ; //bins in ptpi2 from 8 to 10 GeV
163         const Int_t nbinzvtx  = 30 ; //bins in z vertex
164         const Int_t nbincent = 18;  //bins in centrality
165         const Int_t nbincent_0_10 = 4;  //bins in centrality between 0 and 10
166         const Int_t nbincent_10_60 = 10;  //bins in centrality between 10 and 60
167         const Int_t nbincent_60_100 = 4;  //bins in centrality between 60 and 100
168         const Int_t nbinfake = 3;  //bins in fake
169         const Int_t nbindist12 = 10; //bins dist12
170         const Int_t nbindist23 = 10; //bins dist23
171         const Int_t nbinsigmaVtx = 10; //bin sigmaVtx   
172         const Int_t nbinsumd02 = 10; //bin sumD0^2      
173         const Int_t nbinpointingXY = 50;  //bins in cosPointingAngleXY
174         const Int_t nbinnormDecayLXY = 20;  //bins in NormDecayLengthXY
175         const Int_t nbinmult = 48;  //bins in multiplicity (total number)
176         const Int_t nbinmult_0_20 = 20; //bins in multiplicity between 0 and 20
177         const Int_t nbinmult_20_50 = 15; //bins in multiplicity between 20 and 50
178         const Int_t nbinmult_50_102 = 13; //bins in multiplicity between 50 and 102
179         
180         //the sensitive variables, their indices
181         const UInt_t ipT = 0;
182         const UInt_t iy  = 1;
183         const UInt_t iphi  = 2;
184         const UInt_t icT  = 3;
185         const UInt_t ipointing  = 4;
186         const UInt_t ipTpi  = 5;
187         const UInt_t ipTk  = 6;
188         const UInt_t ipTpi2  = 7;
189         const UInt_t izvtx  = 8;
190         const UInt_t icent = 9;
191         const UInt_t ifake = 10;
192         const UInt_t idist12 = 11;
193         const UInt_t idist23 = 12;
194         const UInt_t isigmaVtx = 13;
195         const UInt_t isumd02 = 14;
196         const UInt_t ipointingXY = 15;
197         const UInt_t inormDecayLXY = 16;
198         const UInt_t imult = 17;
199
200         const Int_t nvarTot   = 18 ; //number of variables on the grid:pt, y, cosThetaStar, pTpi, pTk, cT, dca, d0pi, d0K, d0xd0, cosPointingAngle, phi, zvtx, centrality, fake, cosPointingAngleXY, normDecayLengthXY, multiplicity
201
202         //arrays for the number of bins in each dimension
203         Int_t iBin[nvarTot];
204         //iBin[ipT]=nbinpt_0_4+nbinpt_4_8+nbinpt_8_10;
205         iBin[ipT]=nbinpt;
206         iBin[iy]=nbiny;
207         iBin[iphi]=nbinphi;
208         //      iBin[icT]=nbincT_0_4+nbincT_4_8+nbincT_8_10;
209         //iBin[4]=nbinpointing_0_4+nbinpointing_4_8+nbinpointing_8_10;
210         iBin[icT]=nbincT;
211         iBin[ipointing]=nbinpointing;
212         iBin[ipTpi]=nbinpt;
213         iBin[ipTk]=nbinpt;
214         iBin[ipTpi2]=nbinpt;
215         iBin[izvtx]=nbinzvtx;
216         iBin[icent]=nbincent;
217         iBin[ifake]=nbinfake;
218         iBin[idist12]=nbindist12;
219         iBin[idist23]=nbindist23;
220         iBin[isigmaVtx]=nbinsigmaVtx;
221         iBin[isumd02]=nbinsumd02;
222         iBin[ipointingXY]=nbinpointingXY;
223         iBin[inormDecayLXY]=nbinnormDecayLXY;
224         iBin[imult]=nbinmult;
225         
226         //arrays for lower bounds :
227         Double_t *binLimpT=new Double_t[iBin[ipT]+1];
228         Double_t *binLimy=new Double_t[iBin[iy]+1];
229         Double_t *binLimphi=new Double_t[iBin[iphi]+1];
230         Double_t *binLimcT=new Double_t[iBin[icT]+1];
231         Double_t *binLimpointing=new Double_t[iBin[ipointing]+1];
232         Double_t *binLimpTpi=new Double_t[iBin[ipTpi]+1];
233         Double_t *binLimpTk=new Double_t[iBin[ipTk]+1];
234         Double_t *binLimpTpi2=new Double_t[iBin[ipTpi2]+1];
235         Double_t *binLimzvtx=new Double_t[iBin[izvtx]+1];
236         Double_t *binLimcent=new Double_t[iBin[icent]+1];
237         Double_t *binLimfake=new Double_t[iBin[ifake]+1];
238         Double_t *binLimdist12=new Double_t[iBin[idist12]+1];
239         Double_t *binLimdist23=new Double_t[iBin[idist23]+1];
240         Double_t *binLimsigmaVtx=new Double_t[iBin[isigmaVtx]+1];
241         Double_t *binLimsumd02=new Double_t[iBin[isumd02]+1];   
242         Double_t *binLimpointingXY=new Double_t[iBin[ipointingXY]+1];
243         Double_t *binLimnormDecayLXY=new Double_t[iBin[inormDecayLXY]+1];
244         Double_t *binLimmult=new Double_t[iBin[imult]+1];
245         
246         // checking limits
247         /*
248           if (ptmax_0_4 != ptmin_4_8) {
249           Error("AliCFHeavyFlavourTaskMultiVarMultiStep","max lim 1st range != min lim 2nd range, please check!");
250           }
251           if (ptmax_4_8 != ptmin_8_10) {
252           Error("AliCFHeavyFlavourTaskMultiVarMultiStep","max lim 2nd range != min lim 3rd range, please check!");
253           }
254         */
255         // values for bin lower bounds
256         // pt
257         Float_t* floatbinLimpT = cutsLctopKpi->GetPtBinLimits();
258         for (Int_t ibinpT = 0 ; ibinpT<iBin[ipT]+1; ibinpT++){
259                 binLimpT[ibinpT] = (Double_t)floatbinLimpT[ibinpT];
260                 binLimpTpi[ibinpT] = (Double_t)floatbinLimpT[ibinpT];
261                 binLimpTk[ibinpT] = (Double_t)floatbinLimpT[ibinpT];
262                 binLimpTpi2[ibinpT] = (Double_t)floatbinLimpT[ibinpT];
263         }
264         for(Int_t i=0; i<=nbinpt; i++) printf("binLimpT[%d]=%f\n",i,binLimpT[i]);  
265         
266         /*
267           for(Int_t i=0; i<=nbinpt_0_4; i++) binLimpT[i]=(Double_t)ptmin_0_4 + (ptmax_0_4-ptmin_0_4)/nbinpt_0_4*(Double_t)i ; 
268           if (binLimpT[nbinpt_0_4] != ptmin_4_8)  {
269           Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for pt - 1st range - differs from expected!\n");
270           }
271           for(Int_t i=0; i<=nbinpt_4_8; i++) binLimpT[i+nbinpt_0_4]=(Double_t)ptmin_4_8 + (ptmax_4_8-ptmin_4_8)/nbinpt_4_8*(Double_t)i ; 
272           if (binLimpT[nbinpt_0_4+nbinpt_4_8] != ptmin_8_10)  {
273           Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for pt - 2nd range - differs from expected!\n");
274           }
275           for(Int_t i=0; i<=nbinpt_8_10; i++) binLimpT[i+nbinpt_0_4+nbinpt_4_8]=(Double_t)ptmin_8_10 + (ptmax_8_10-ptmin_8_10)/nbinpt_8_10*(Double_t)i ; 
276         */
277         
278         // y
279         for(Int_t i=0; i<=nbiny; i++) binLimy[i]=(Double_t)ymin  + (ymax-ymin)  /nbiny*(Double_t)i ;
280         
281         // Phi
282         for(Int_t i=0; i<=nbinphi; i++) binLimphi[i]=(Double_t)phimin  + (phimax-phimin)  /nbinphi*(Double_t)i ;
283         
284         // cT
285         for(Int_t i=0; i<=nbincT; i++) binLimcT[i]=(Double_t)cTmin  + (cTmax-cTmin)  /nbincT*(Double_t)i ;
286         
287         // cosPointingAngle
288         for(Int_t i=0; i<=nbinpointing; i++) binLimpointing[i]=(Double_t)cosmin  + (cosmax-cosmin)  /nbinpointing*(Double_t)i ;
289
290         /*
291         // ptPi
292         for(Int_t i=0; i<=nbincT_0_4; i++) binLimcT[i]=(Double_t)ptmin_0_4 + (ptmax_0_4-ptmin_0_4)/nbincT_0_4*(Double_t)i ; 
293         if (binLimcT[nbincT_0_4] != ptmin_4_8)  {
294         Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for ptPi - 1st range - differs from expected!");
295         }
296         for(Int_t i=0; i<=nbincT_4_8; i++) binLimcT[i+nbincT_0_4]=(Double_t)ptmin_4_8 + (ptmax_4_8-ptmin_4_8)/nbincT_4_8*(Double_t)i ; 
297         if (binLimcT[nbincT_0_4+nbincT_4_8] != ptmin_8_10)  {
298         Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for ptPi - 2nd range - differs from expected!\n");
299         }
300         for(Int_t i=0; i<=nbincT_8_10; i++) binLimcT[i+nbincT_0_4+nbincT_4_8]=(Double_t)ptmin_8_10 + (ptmax_8_10-ptmin_8_10)/nbincT_8_10*(Double_t)i ; 
301         
302         // ptKa
303         for(Int_t i=0; i<=nbinpointing_0_4; i++) binLimpointing[i]=(Double_t)ptmin_0_4 + (ptmax_0_4-ptmin_0_4)/nbinpointing_0_4*(Double_t)i ; 
304         if (binLimpointing[nbinpointing_0_4] != ptmin_4_8)  {
305         Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for ptKa - 1st range - differs from expected!");
306         }
307         for(Int_t i=0; i<=nbinpointing_4_8; i++) binLimpointing[i+nbinpointing_0_4]=(Double_t)ptmin_4_8 + (ptmax_4_8-ptmin_4_8)/nbinpointing_4_8*(Double_t)i ; 
308         if (binLimpointing[nbinpointing_0_4+nbinpointing_4_8] != ptmin_8_10)  {
309         Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for ptKa - 2nd range - differs from expected!\n");
310         }
311         for(Int_t i=0; i<=nbinpointing_8_10; i++) binLimpointing[i+nbinpointing_0_4+nbinpointing_4_8]=(Double_t)ptmin_8_10 + (ptmax_8_10-ptmin_8_10)/nbinpointing_8_10*(Double_t)i ; 
312         */
313         
314         // z Primary Vertex
315         for(Int_t i=0; i<=nbinzvtx; i++) {
316                 binLimzvtx[i]=(Double_t)zvtxmin  + (zvtxmax-zvtxmin)  /nbinzvtx*(Double_t)i ;
317         }
318         
319         // centrality
320         for(Int_t i=0; i<=nbincent_0_10; i++) binLimcent[i]=(Double_t)centmin_0_10 + (centmax_0_10-centmin_0_10)/nbincent_0_10*(Double_t)i ; 
321         if (binLimcent[nbincent_0_10] != centmin_10_60)  {
322                 Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for cent - 1st range - differs from expected!\n");
323         }
324         for(Int_t i=0; i<=nbincent_10_60; i++) binLimcent[i+nbincent_0_10]=(Double_t)centmin_10_60 + (centmax_10_60-centmin_10_60)/nbincent_10_60*(Double_t)i ;
325         if (binLimcent[nbincent_0_10+nbincent_10_60] != centmin_60_100)  {
326                 Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for cent - 2st range - differs from expected!\n");
327         }
328         for(Int_t i=0; i<=nbincent_60_100; i++) binLimcent[i+nbincent_10_60]=(Double_t)centmin_60_100 + (centmax_60_100-centmin_60_100)/nbincent_60_100*(Double_t)i ;
329         
330         // fake
331         for(Int_t i=0; i<=nbinfake; i++) {
332           binLimfake[i]=(Double_t)fakemin  + (fakemax-fakemin)/nbinfake * (Double_t)i;
333         }
334
335         //dist12
336         for(Int_t i=0; i<=nbindist12; i++) {
337           binLimdist12[i]=(Double_t)distTwoPartmin  + (distTwoPartmax-distTwoPartmin)/nbindist12 * (Double_t)i;
338         }
339  
340         //dist23
341         for(Int_t i=0; i<=nbindist23; i++) {
342           binLimdist23[i]=(Double_t)distTwoPartmin  + (distTwoPartmax-distTwoPartmin)/nbindist23 * (Double_t)i;
343         }
344  
345         //dispersion Vtx
346         for(Int_t i=0; i<=nbinsigmaVtx; i++) {
347           binLimsigmaVtx[i]=(Double_t)dispVtxmin  + (dispVtxmax-dispVtxmin)/nbinsigmaVtx * (Double_t)i;
348         }
349
350         //sumd0^2
351         for(Int_t i=0; i<=nbinsumd02; i++) {
352           binLimsumd02[i]=(Double_t)sumd02min  + (sumd02max-sumd02min)/nbinsumd02 * (Double_t)i;
353         }
354
355         // cosPointingAngleXY
356         for(Int_t i=0; i<=nbinpointingXY; i++) binLimpointingXY[i]=(Double_t)cosminXY  + (cosmaxXY-cosminXY)  /nbinpointingXY*(Double_t)i ;
357
358         // normDecayLXY
359         for(Int_t i=0; i<=nbinnormDecayLXY; i++) binLimnormDecayLXY[i]=(Double_t)normDecLXYmin  + (normDecLXYmax-normDecLXYmin)  /nbinnormDecayLXY*(Double_t)i ;
360
361         // multiplicity
362         for(Int_t i=0; i<=nbinmult_0_20; i++) binLimmult[i]=(Double_t)multmin_0_20 + (multmax_0_20-multmin_0_20)/nbinmult_0_20*(Double_t)i ; 
363         if (binLimmult[nbinmult_0_20] != multmin_20_50)  {
364                 Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for mult - 1st range - differs from expected!\n");
365         }
366         for(Int_t i=0; i<=nbinmult_20_50; i++) binLimmult[i+nbinmult_0_20]=(Double_t)multmin_20_50 + (multmax_20_50-multmin_20_50)/nbinmult_20_50*(Double_t)i ; 
367         if (binLimmult[nbinmult_0_20+nbinmult_20_50] != multmin_50_102)  {
368                 Error("AliCFHeavyFlavourTaskMultiVarMultiStep","Calculated bin lim for mult - 2nd range - differs from expected!\n");
369         }
370         for(Int_t i=0; i<=nbinmult_50_102; i++) binLimmult[i+nbinmult_0_20+nbinmult_20_50]=(Double_t)multmin_50_102 + (multmax_50_102-multmin_50_102)/nbinmult_50_102*(Double_t)i ; 
371         
372         //one "container" for MC
373         TString nameContainer="";
374         if(!isKeepDfromB) {
375                 nameContainer="CFHFccontainer0_3Prong_CommonFramework";
376         }
377         else  if(isKeepDfromBOnly){
378                 nameContainer="CFHFccontainer0DfromB_3Prong_CommonFramework";
379         }
380         else  {
381                 nameContainer="CFHFccontainer0allD_3Prong_CommonFramework";          
382         }
383         
384         AliCFContainer* container;
385         if (configuration == AliCFTaskVertexingHF::kSnail){
386                 container = new AliCFContainer(nameContainer,"container for tracks",nstep,nvarTot,iBin);
387                 //setting the bin limits
388                 printf("pt\n");
389                 container -> SetBinLimits(ipT,binLimpT);
390                 printf("y\n");
391                 container -> SetBinLimits(iy,binLimy);
392                 printf("Phi\n");
393                 container -> SetBinLimits(iphi,binLimphi);
394                 printf("cT\n");
395                 container -> SetBinLimits(icT,binLimcT);
396                 printf("pointing angle\n");
397                 container -> SetBinLimits(ipointing,binLimpointing);
398                 printf("ptpi\n");
399                 container -> SetBinLimits(ipTpi,binLimpTpi);
400                 printf("ptK\n");
401                 container -> SetBinLimits(ipTk,binLimpTk);
402                 printf("ptpi2\n");
403                 container -> SetBinLimits(ipTpi2,binLimpTpi2);
404                 printf("zvtx \n");
405                 container -> SetBinLimits(izvtx,binLimzvtx);
406                 printf("cent\n");
407                 container -> SetBinLimits(icent,binLimcent);
408                 printf("fake\n");
409                 container -> SetBinLimits(ifake,binLimfake);
410                 printf("dist12\n");
411                 container -> SetBinLimits(idist12,binLimdist12);
412                 printf("dist23\n");
413                 container -> SetBinLimits(idist23,binLimdist23);
414                 printf("dispVtx\n");
415                 container -> SetBinLimits(isigmaVtx,binLimsigmaVtx);
416                 printf("sumd0^2\n");
417                 container -> SetBinLimits(isumd02,binLimsumd02);
418                 printf("pointingXY\n");
419                 container -> SetBinLimits(ipointingXY,binLimpointingXY);
420                 printf("normDecayLXY\n");
421                 container -> SetBinLimits(inormDecayLXY,binLimnormDecayLXY);
422                 printf("multiplicity\n");
423                 container -> SetBinLimits(imult,binLimmult);
424                 
425                 container -> SetVarTitle(ipT,"pt");
426                 container -> SetVarTitle(iy,"y");
427                 container -> SetVarTitle(iphi, "phi");
428                 container -> SetVarTitle(icT, "ct");
429                 container -> SetVarTitle(ipointing, "pionting");        
430                 container -> SetVarTitle(ipTpi, "ptpi");
431                 container -> SetVarTitle(ipTk, "ptK");
432                 container -> SetVarTitle(ipTpi2, "ptpi2");
433                 container -> SetVarTitle(izvtx, "zvtx");
434                 container -> SetVarTitle(icent, "centrality");
435                 container -> SetVarTitle(ifake, "fake");
436                 container -> SetVarTitle(idist12, "dist12toVtx");
437                 container -> SetVarTitle(idist23, "dist23toVtx");
438                 container -> SetVarTitle(isigmaVtx, "dispertionToSecVtx");
439                 container -> SetVarTitle(isumd02, "sumd0^2");
440                 container -> SetVarTitle(ipointingXY, "piointingXY");
441                 container -> SetVarTitle(inormDecayLXY, "normDecayLXY");
442                 container -> SetVarTitle(imult, "multiplicity");
443         }
444         else if (configuration == AliCFTaskVertexingHF::kCheetah){
445                 //arrays for the number of bins in each dimension
446                 const Int_t nvar = 8;
447
448                 const UInt_t ipTFast = 0;
449                 const UInt_t iyFast = 1;
450                 const UInt_t icTFast = 2;
451                 const UInt_t iphiFast = 3;
452                 const UInt_t izvtxFast = 4;
453                 const UInt_t icentFast = 5;
454                 const UInt_t ifakeFast = 6;
455                 const UInt_t imultFast = 7;
456
457                 Int_t iBinFast[nvar];
458                 iBinFast[ipTFast] = iBin[ipT];
459                 iBinFast[iyFast] = iBin[iy];
460                 iBinFast[icTFast] = iBin[icT];
461                 iBinFast[iphiFast] = iBin[iphi];
462                 iBinFast[izvtxFast] = iBin[izvtx];
463                 iBinFast[icentFast] = iBin[icent];
464                 iBinFast[ifakeFast] = iBin[ifake];
465                 iBinFast[imultFast] = iBin[imult];
466
467                 container = new AliCFContainer(nameContainer,"container for tracks",nstep,nvar,iBinFast);
468                 printf("pt\n");
469                 container -> SetBinLimits(ipTFast,binLimpT);
470                 printf("y\n");
471                 container -> SetBinLimits(iyFast,binLimy);
472                 printf("ct\n");
473                 container -> SetBinLimits(icTFast,binLimcT);
474                 printf("phi\n");
475                 container -> SetBinLimits(iphiFast,binLimphi);
476                 printf("zvtx\n");
477                 container -> SetBinLimits(izvtxFast,binLimzvtx);
478                 printf("centrality\n");
479                 container -> SetBinLimits(icentFast,binLimcent);
480                 printf("fake\n");
481                 container -> SetBinLimits(ifakeFast,binLimfake);
482                 printf("multiplicity\n");
483                 container -> SetBinLimits(imultFast,binLimmult);
484
485                 container -> SetVarTitle(ipTFast,"pt");
486                 container -> SetVarTitle(iyFast,"y");
487                 container -> SetVarTitle(icTFast, "ct");
488                 container -> SetVarTitle(iphiFast, "phi");
489                 container -> SetVarTitle(izvtxFast, "zvtx");
490                 container -> SetVarTitle(icentFast, "centrality");
491                 container -> SetVarTitle(ifakeFast, "fake");
492                 container -> SetVarTitle(imultFast, "multiplicity");
493         }
494
495         return container;
496
497         container -> SetStepTitle(0, "MCLimAcc");
498         container -> SetStepTitle(1, "MC");
499         container -> SetStepTitle(2, "MCAcc");
500         container -> SetStepTitle(3, "RecoVertex");
501         container -> SetStepTitle(4, "RecoRefit");
502         container -> SetStepTitle(5, "Reco");
503         container -> SetStepTitle(6, "RecoAcc");
504         container -> SetStepTitle(7, "RecoITSCluster");
505         container -> SetStepTitle(8, "RecoCuts");
506         container -> SetStepTitle(9, "RecoPID");
507
508
509         //CREATE THE  CUTS -----------------------------------------------
510         
511         // Gen-Level kinematic cuts
512         AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
513         
514         //Particle-Level cuts:  
515         AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
516         Bool_t useAbsolute = kTRUE;
517         if (isSign != 2){
518                 useAbsolute = kFALSE;
519         }
520         mcGenCuts->SetRequirePdgCode(pdgCode, useAbsolute);  // kTRUE set in order to include antiparticle
521         mcGenCuts->SetAODMC(1); //special flag for reading MC in AOD tree (important)
522         
523         // Acceptance cuts:
524         AliCFAcceptanceCuts* accCuts = new AliCFAcceptanceCuts("accCuts", "Acceptance cuts");
525         AliCFTrackKineCuts *kineAccCuts = new AliCFTrackKineCuts("kineAccCuts","Kine-Acceptance cuts");
526         kineAccCuts->SetPtRange(ptmin,ptmax);
527         kineAccCuts->SetEtaRange(etamin,etamax);
528
529         // Rec-Level kinematic cuts
530         AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
531         
532         AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
533         
534         AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
535         
536         printf("CREATE MC KINE CUTS\n");
537         TObjArray* mcList = new TObjArray(0) ;
538         mcList->AddLast(mcKineCuts);
539         mcList->AddLast(mcGenCuts);
540         
541         printf("CREATE ACCEPTANCE CUTS\n");
542         TObjArray* accList = new TObjArray(0) ;
543         accList->AddLast(kineAccCuts);
544
545         printf("CREATE RECONSTRUCTION CUTS\n");
546         TObjArray* recList = new TObjArray(0) ;   // not used!! 
547         recList->AddLast(recKineCuts);
548         recList->AddLast(recQualityCuts);
549         recList->AddLast(recIsPrimaryCuts);
550         
551         TObjArray* emptyList = new TObjArray(0);
552
553         //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
554         printf("CREATE INTERFACE AND CUTS\n");
555         AliCFManager* man = new AliCFManager() ;
556         man->SetParticleContainer(container);
557         man->SetParticleCutsList(0 , mcList); // MC, Limited Acceptance
558         man->SetParticleCutsList(1 , mcList); // MC
559         man->SetParticleCutsList(2 , accList); // Acceptance 
560         man->SetParticleCutsList(3 , emptyList); // Vertex 
561         man->SetParticleCutsList(4 , emptyList); // Refit 
562         man->SetParticleCutsList(5 , emptyList); // AOD
563         man->SetParticleCutsList(6 , emptyList); // AOD in Acceptance
564         man->SetParticleCutsList(7 , emptyList); // AOD with required n. of ITS clusters
565         man->SetParticleCutsList(8 , emptyList); // AOD Reco (PPR cuts implemented in Task)
566         man->SetParticleCutsList(9 , emptyList); // AOD Reco PID
567         
568         // Get the pointer to the existing analysis manager via the static access method.
569         //==============================================================================
570         AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
571         if (!mgr) {
572           ::Error("AddTaskCompareHF", "No analysis manager to connect to.");
573           return NULL;
574         }   
575         //CREATE THE TASK
576         printf("CREATE TASK\n");
577
578         // create the task
579         AliCFTaskVertexingHF *task = new AliCFTaskVertexingHF("AliCFTaskVertexingHF",cutsLctopKpi);
580         task->SetFillFromGenerated(kFALSE);
581         task->SetDecayChannel(32);
582         task->SetUseWeight(kFALSE);
583         task->SetCFManager(man); //here is set the CF manager
584         task->SetSign(isSign);
585         task->SetCentralitySelection(kFALSE);
586         task->SetFakeSelection(0);
587         task->SetRejectCandidateIfNotFromQuark(kTRUE); // put to false if you want to keep HIJING D0!!
588         task->SetUseMCVertex(kFALSE); // put to true if you want to do studies on pp
589         if (isKeepDfromB && !isKeepDfromBOnly) task->SetDselection(2);
590         if (isKeepDfromB && isKeepDfromBOnly) task->SetDselection(1);           
591
592         TF1* funcWeight = 0x0;
593         if (task->GetUseWeight()) {
594                 funcWeight = (TF1*)cutFile->Get("funcWeight");
595                 if (funcWeight == 0x0){
596                         Printf("FONLL Weights will be used");
597                 }
598                 else {
599                         task->SetWeightFunction(funcWeight);
600                         Printf("User-defined Weights will be used. The function being:");
601                         task->GetWeightFunction()->Print();
602                 }
603         }
604
605         Printf("***************** CONTAINER SETTINGS *****************");
606         Printf("decay channel = %d",(Int_t)task->GetDecayChannel());
607         Printf("FillFromGenerated = %d",(Int_t)task->GetFillFromGenerated());
608         Printf("Dselection = %d",(Int_t)task->GetDselection());
609         Printf("UseWeight = %d",(Int_t)task->GetUseWeight());
610         if (task->GetUseWeight()) {
611                 funcWeight = (TF1*)cutFile->Get("funcWeight");
612                 if (funcWeight == 0x0){
613                         Printf("FONLL Weights will be used");
614                 }
615                 else {
616                         task->SetWeightFunction(funcWeight);
617                         Printf("User-defined Weights will be used. The function being:");
618                         task->GetWeightFunction()->Print();
619                 }
620         }
621         Printf("Sign = %d",(Int_t)task->GetSign());
622         Printf("Centrality selection = %d",(Int_t)task->GetCentralitySelection());
623         Printf("Fake selection = %d",(Int_t)task->GetFakeSelection());
624         Printf("RejectCandidateIfNotFromQuark selection = %d",(Int_t)task->GetRejectCandidateIfNotFromQuark());
625         Printf("UseMCVertex selection = %d",(Int_t)task->GetUseMCVertex());
626         Printf("***************END CONTAINER SETTINGS *****************\n");
627
628         //-----------------------------------------------------------//
629         //   create correlation matrix for unfolding - only eta-pt   //
630         //-----------------------------------------------------------//
631
632         Bool_t AcceptanceUnf = kTRUE; // unfold at acceptance level, otherwise PPR
633
634         Int_t thnDim[4];
635         
636         //first half  : reconstructed 
637         //second half : MC
638
639         thnDim[0] = iBin[ipT];
640         thnDim[2] = iBin[ipT];
641         thnDim[1] = iBin[iy];
642         thnDim[3] = iBin[iy];
643
644         TString nameCorr="";
645         if(!isKeepDfromB) {
646                 nameCorr="CFHFcorr0_3Prong_CommonFramework";
647         }
648         else  if(isKeepDfromBOnly){
649                 nameCorr= "CFHFcorr0KeepDfromBOnly_3Prong_CommonFramework";
650         }
651         else  {
652                 nameCorr="CFHFcorr0allD_3Prong_CommonFramework";                
653         }
654
655         THnSparseD* correlation = new THnSparseD(nameCorr,"THnSparse with correlations",4,thnDim);
656         Double_t** binEdges = new Double_t[2];
657
658         // set bin limits
659
660         binEdges[0]= binLimpT;
661         binEdges[1]= binLimy;
662
663         correlation->SetBinEdges(0,binEdges[0]);
664         correlation->SetBinEdges(2,binEdges[0]);
665
666         correlation->SetBinEdges(1,binEdges[1]);
667         correlation->SetBinEdges(3,binEdges[1]);
668
669         correlation->Sumw2();
670   
671         // correlation matrix ready
672         //------------------------------------------------//
673
674         task->SetCorrelationMatrix(correlation); // correlation matrix for unfolding
675         
676         // Create and connect containers for input/output
677         
678         // ------ input data ------
679         AliAnalysisDataContainer *cinput0  = mgr->GetCommonInputContainer();
680         
681         // ----- output data -----
682         
683         TString outputfile = AliAnalysisManager::GetCommonFileName();
684         TString output1name="", output2name="", output3name="", output4name="";;
685         output2name=nameContainer;
686         output3name=nameCorr;
687         if(!isKeepDfromB) {
688                 outputfile += ":PWG3_D2H_CFtaskLctopKpi_CommonFramework";
689                 output1name="CFHFchist0_3Prong_CommonFramework";
690         }
691         else  if(isKeepDfromBOnly){
692                 outputfile += ":PWG3_D2H_CFtaskLctopKpiKeepDfromBOnly_CommonFramework";
693                 output1name="CFHFchist0DfromB_3Prong_CommonFramework";
694         }
695         else{
696                 outputfile += ":PWG3_D2H_CFtaskLctopKpiKeepDfromB_CommonFramework";
697                 output1name="CFHFchist0allD_3Prong_CommonFramework";
698         }
699
700         output4name= "Cuts_3Prong_CommonFramework";
701
702         //now comes user's output objects :
703         // output TH1I for event counting
704         AliAnalysisDataContainer *coutput1 = mgr->CreateContainer(output1name, TH1I::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
705         // output Correction Framework Container (for acceptance & efficiency calculations)
706         AliAnalysisDataContainer *coutput2 = mgr->CreateContainer(output2name, AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
707         // Unfolding - correlation matrix
708         AliAnalysisDataContainer *coutput3 = mgr->CreateContainer(output3name, THnSparseD::Class(),AliAnalysisManager::kOutputContainer,outputfile.Data());
709         AliAnalysisDataContainer *coutput4 = mgr->CreateContainer(output4name, AliRDHFCuts::Class(),AliAnalysisManager::kOutputContainer, outputfile.Data());
710
711         mgr->AddTask(task);
712         
713         mgr->ConnectInput(task,0,mgr->GetCommonInputContainer());
714         mgr->ConnectOutput(task,1,coutput1);
715         mgr->ConnectOutput(task,2,coutput2);
716         mgr->ConnectOutput(task,3,coutput3);
717         mgr->ConnectOutput(task,4,coutput4);
718
719         return task;
720         }
721