Added option to set event selection
[u/mrichter/AliRoot.git] / PWGCF / FLOW / macros / AddTaskFlowCentralityPIDSP.C
1 /////////////////////////////////////////////////////////////////////////////////////////////\r
2 //\r
3 // AddTask* macro for flow analysis\r
4 // Creates a Flow Event task and adds it to the analysis manager.\r
5 // Sets the cuts using the correction framework (CORRFW) classes.\r
6 // Also creates Flow Analysis tasks and connects them to the output of the flow event task.\r
7 //\r
8 /////////////////////////////////////////////////////////////////////////////////////////////\r
9 \r
10 void AddTaskFlowCentralityPIDSP(Int_t centralitysel,\r
11                                Float_t centrMin=10.,\r
12                                Float_t centrMax=20.,\r
13                                TString fileNameBase="output",\r
14                                Bool_t isPID = kTRUE,\r
15                                AliPID::EParticleType particleType=AliPID::kPion,\r
16                                AliFlowTrackCuts::PIDsource sourcePID = AliFlowTrackCuts::kTOFbayesian,\r
17                                Int_t charge=0,\r
18                                Int_t harmonic=2,\r
19                                Bool_t doQA=kFALSE,\r
20                                Float_t etamin=-0.8,\r
21                                Float_t etamax=0.8,       \r
22                                TString uniqueStr="" )\r
23 {\r
24   // Define the range for eta subevents (for SP method)\r
25   Double_t minA = -5;\r
26   Double_t maxA = -1.5;\r
27   Double_t minB = 1.5;\r
28   Double_t maxB = 5;\r
29 \r
30   // AFTERBURNER\r
31   Bool_t useAfterBurner=kFALSE;\r
32   Double_t v1=0.0;\r
33   Double_t v2=0.0;\r
34   Double_t v3=0.0;\r
35   Double_t v4=0.0;\r
36   Int_t numberOfTrackClones=0; //non-flow\r
37 \r
38   // Define a range of the detector to exclude\r
39   Bool_t ExcludeRegion = kFALSE;\r
40   Double_t excludeEtaMin = -0.;\r
41   Double_t excludeEtaMax = 0.;\r
42   Double_t excludePhiMin = 0.;\r
43   Double_t excludePhiMax = 0.;\r
44 \r
45   // use physics selection class\r
46   Bool_t  UsePhysicsSelection = kTRUE;\r
47 \r
48   // QA\r
49   Bool_t runQAtask=kFALSE;\r
50   Bool_t FillQAntuple=kFALSE;\r
51   Bool_t DoQAcorrelations=kFALSE;\r
52 \r
53   // RUN SETTINGS\r
54   // Flow analysis method can be:(set to kTRUE or kFALSE)\r
55   Bool_t SP       = kTRUE;  // scalar product method (similar to eventplane method)\r
56   Bool_t QC       = kFALSE;  // cumulants using Q vectors\r
57   \r
58   //these are OBSOLETE, use at own peril\r
59   Bool_t GFC      = kFALSE;  // cumulants based on generating function\r
60   Bool_t MCEP     = kFALSE;  // correlation with Monte Carlo reaction plane\r
61   Bool_t FQD      = kFALSE;  // fit of the distribution of the Q vector (only integrated v)\r
62   Bool_t LYZ1SUM  = kFALSE;  // Lee Yang Zeroes using sum generating function (integrated v)\r
63   Bool_t LYZ1PROD = kFALSE;  // Lee Yang Zeroes using product generating function (integrated v)\r
64   Bool_t LYZ2SUM  = kFALSE; // Lee Yang Zeroes using sum generating function (second pass differential v)\r
65   Bool_t LYZ2PROD = kFALSE; // Lee Yang Zeroes using product generating function (second pass differential v)\r
66   Bool_t LYZEP    = kFALSE; // Lee Yang Zeroes Event plane using sum generating function (gives eventplane + weight)\r
67   Bool_t MH       = kFALSE;  // azimuthal correlators in mixed harmonics  \r
68   Bool_t NL       = kFALSE;  // nested loops (for instance distribution of phi1-phi2 for all distinct pairs)\r
69 \r
70   Bool_t METHODS[] = {SP,LYZ1SUM,LYZ1PROD,LYZ2SUM,LYZ2PROD,LYZEP,GFC,QC,FQD,MCEP,MH,NL};\r
71 \r
72   // Boolean to use/not use weights for the Q vector\r
73   Bool_t WEIGHTS[] = {kFALSE,kFALSE,kFALSE}; //Phi, v'(pt), v'(eta)\r
74 \r
75   // SETTING THE CUTS\r
76 \r
77   //---------Data selection----------\r
78   //kMC, kGlobal, kESD_TPConly, kESD_SPDtracklet\r
79   AliFlowTrackCuts::trackParameterType rptype = AliFlowTrackCuts::kTPCstandalone;\r
80   AliFlowTrackCuts::trackParameterType poitype = AliFlowTrackCuts::kTPCstandalone;\r
81 \r
82   //---------Parameter mixing--------\r
83   //kPure - no mixing, kTrackWithMCkine, kTrackWithMCPID, kTrackWithMCpt\r
84   AliFlowTrackCuts::trackParameterMix rpmix = AliFlowTrackCuts::kPure;\r
85   AliFlowTrackCuts::trackParameterMix poimix = AliFlowTrackCuts::kPure;\r
86 \r
87 \r
88   const char* rptypestr = AliFlowTrackCuts::GetParamTypeName(rptype);\r
89   const char* poitypestr = AliFlowTrackCuts::GetParamTypeName(poitype);\r
90 \r
91   //===========================================================================\r
92   // EVENTS CUTS:\r
93  AliFlowEventCuts* cutsEvent = new AliFlowEventCuts("event cuts");\r
94   cutsEvent->SetCentralityPercentileRange(centrMin,centrMax);\r
95   cutsEvent->SetCentralityPercentileMethod(AliFlowEventCuts::kV0);\r
96 //  cutsEvent->SetRefMultMethod(AliFlowEventCuts::kV0);\r
97   //cutsEvent->SetCentralityPercentileMethod(AliFlowEventCuts::kSPD1tracklets);\r
98 //  cutsEvent->SetNContributorsRange(2);\r
99   cutsEvent->SetPrimaryVertexZrange(-10.,10.);\r
100   cutsEvent->SetQA(doQA);\r
101   cutsEvent->SetCutTPCmultiplicityOutliers();\r
102   \r
103   // RP TRACK CUTS:\r
104   AliFlowTrackCuts* cutsRP2 = AliFlowTrackCuts::GetStandardVZEROOnlyTrackCuts();\r
105   AliFlowTrackCuts* cutsRP = new AliFlowTrackCuts("TPConlyRP");\r
106   cutsRP->SetParamType(rptype);\r
107   cutsRP->SetParamMix(rpmix);\r
108   cutsRP->SetPtRange(0.2,5.);\r
109   cutsRP->SetEtaRange(etamin,etamax);\r
110   cutsRP->SetMinNClustersTPC(70);\r
111 //  cutsRP->SetMinChi2PerClusterTPC(0.1);\r
112 //  cutsRP->SetMaxChi2PerClusterTPC(4.0);\r
113   cutsRP->SetMaxDCAToVertexXY(3.0);\r
114   cutsRP->SetMaxDCAToVertexZ(3.0);\r
115   cutsRP->SetAcceptKinkDaughters(kFALSE);\r
116   cutsRP->SetMinimalTPCdedx(10.);\r
117   cutsRP->SetQA(doQA);\r
118 \r
119   // POI TRACK CUTS:\r
120   AliFlowTrackCuts* cutsPOI = new AliFlowTrackCuts("TPConlyPOI");\r
121   cutsPOI->SetParamType(poitype);\r
122   cutsPOI->SetParamMix(poimix);\r
123   cutsPOI->SetPtRange(0.0,10.);\r
124   cutsPOI->SetEtaRange(etamin,etamax);\r
125   cutsPOI->SetMinNClustersTPC(70);\r
126 //  cutsPOI->SetMinChi2PerClusterTPC(0.1);\r
127 //  cutsPOI->SetMaxChi2PerClusterTPC(4.0);\r
128 //  cutsPOI->SetRequireITSRefit(kTRUE);\r
129 //  cutsPOI->SetRequireTPCRefit(kTRUE);\r
130 //  cutsPOI->SetMinNClustersITS(2);\r
131   //cutsPOI->SetMaxChi2PerClusterITS(1.e+09);\r
132   cutsPOI->SetMaxDCAToVertexXY(3.0);\r
133   cutsPOI->SetMaxDCAToVertexZ(3.0);\r
134   //cutsPOI->SetDCAToVertex2D(kTRUE);\r
135   //cutsPOI->SetMaxNsigmaToVertex(1.e+10);\r
136   //cutsPOI->SetRequireSigmaToVertex(kFALSE);\r
137   cutsPOI->SetAcceptKinkDaughters(kFALSE);\r
138   if(isPID) cutsPOI->SetPID(particleType, sourcePID);\r
139   if (charge!=0) cutsPOI->SetCharge(charge);\r
140   //cutsPOI->SetAllowTOFmismatch(kFALSE);\r
141   cutsPOI->SetRequireStrictTOFTPCagreement(kTRUE);\r
142   //iexample: francesco's tunig TPC Bethe Bloch for data:\r
143   //cutsPOI->GetESDpid().GetTPCResponse().SetBetheBlochParameters(4.36414e-02,1.75977e+01,1.14385e-08,2.27907e+00,3.36699e+00);\r
144   //cutsPOI->GetESDpid().GetTPCResponse().SetMip(49);\r
145   cutsPOI->SetMinimalTPCdedx(10.);\r
146   cutsPOI->SetAODfilterBit(1);\r
147   cutsPOI->SetQA(doQA);\r
148   cutsPOI->SetPriors((centrMin+centrMax)*0.5); // set priors and PID as a function of the centrality\r
149 \r
150   TString outputSlotName("");\r
151   outputSlotName+=uniqueStr;\r
152   outputSlotName+=Form("V%i ",harmonic);\r
153   outputSlotName+=cutsRP->GetName();\r
154   outputSlotName+=" ";\r
155   outputSlotName+=cutsPOI->GetName();\r
156   outputSlotName+=Form(" %.0f-",centrMin);\r
157   outputSlotName+=Form("%.0f ",centrMax);\r
158   if(isPID){\r
159      outputSlotName+=AliFlowTrackCuts::PIDsourceName(sourcePID);\r
160      outputSlotName+=" ";\r
161      outputSlotName+=AliPID::ParticleName(particleType);\r
162   }\r
163   else{\r
164      outputSlotName+="AllCharged";\r
165   }\r
166   if (charge<0) outputSlotName+="-";\r
167   if (charge>0) outputSlotName+="+";\r
168 \r
169   TString fileName(fileNameBase);\r
170   fileName.Append(".root");\r
171 \r
172   Bool_t useWeights  = WEIGHTS[0] || WEIGHTS[1] || WEIGHTS[2];\r
173   if (useWeights) cout<<"Weights are used"<<endl;\r
174   else cout<<"Weights are not used"<<endl;\r
175   \r
176   // Get the pointer to the existing analysis manager via the static access method.\r
177   //==============================================================================\r
178   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();\r
179   if (!mgr) {\r
180     Error("AddTaskFlowEvent", "No analysis manager to connect to.");\r
181     return NULL;\r
182   }\r
183   \r
184   // Check the analysis type using the event handlers connected to the analysis\r
185   // manager. The availability of MC handler can also be checked here.\r
186   //==============================================================================\r
187   if (!mgr->GetInputEventHandler()) {\r
188     ::Error("AddTaskFlowEvent", "This task requires an input event handler");\r
189     return NULL;\r
190   }  \r
191 \r
192   // Open external input files\r
193   //===========================================================================\r
194   //weights: \r
195   TFile *weightsFile = NULL;\r
196   TList *weightsList = NULL;\r
197 \r
198   if(useWeights) {\r
199     //open the file with the weights:\r
200     weightsFile = TFile::Open("weights.root","READ");\r
201     if(weightsFile) {\r
202       //access the list which holds the histos with weigths:\r
203       weightsList = (TList*)weightsFile->Get("weights");\r
204     }\r
205     else {\r
206       cout<<" WARNING: the file <weights.root> with weights from the previous run was not available."<<endl;\r
207       break;\r
208     } \r
209   }\r
210   \r
211   //LYZ2\r
212   if (LYZ2SUM || LYZ2PROD) {\r
213     //read the outputfile of the first run\r
214     TString outputFileName = "AnalysisResults1.root";\r
215     TString pwd(gSystem->pwd());\r
216     pwd+="/";\r
217     pwd+=outputFileName.Data();\r
218     TFile *outputFile = NULL;\r
219     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {\r
220       cout<<"WARNING: You do not have an output file:"<<endl;\r
221       cout<<"         "<<pwd.Data()<<endl;\r
222       exit(0);\r
223     } else { outputFile = TFile::Open(pwd.Data(),"READ");}\r
224     \r
225     if (LYZ2SUM){  \r
226       // read the output directory from LYZ1SUM \r
227       TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis" ;\r
228       inputFileNameLYZ2SUM += rptypestr;\r
229       cout<<"The input directory is "<<inputFileNameLYZ2SUM.Data()<<endl;\r
230       TFile* fInputFileLYZ2SUM = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2SUM.Data());\r
231       if(!fInputFileLYZ2SUM || fInputFileLYZ2SUM->IsZombie()) { \r
232         cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ; \r
233         break;\r
234       }\r
235       else {\r
236         TList* fInputListLYZ2SUM = (TList*)fInputFileLYZ2SUM->Get("LYZ1SUM");\r
237         if (!fInputListLYZ2SUM) {cout<<"list is NULL pointer!"<<endl;}\r
238       }\r
239       cout<<"LYZ2SUM input file/list read..."<<endl;\r
240     }\r
241 \r
242     if (LYZ2PROD){  \r
243       // read the output directory from LYZ1PROD \r
244       TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis" ;\r
245       inputFileNameLYZ2PROD += rptypestr;\r
246       cout<<"The input directory is "<<inputFileNameLYZ2PROD.Data()<<endl;\r
247       TFile* fInputFileLYZ2PROD = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2PROD.Data());\r
248       if(!fInputFileLYZ2PROD || fInputFileLYZ2PROD->IsZombie()) { \r
249         cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ; \r
250         break;\r
251       }\r
252       else {\r
253         TList* fInputListLYZ2PROD = (TList*)fInputFileLYZ2PROD->Get("LYZ1PROD");\r
254         if (!fInputListLYZ2PROD) {cout<<"list is NULL pointer!"<<endl;}\r
255       }\r
256       cout<<"LYZ2PROD input file/list read..."<<endl;\r
257     }\r
258   }\r
259 \r
260   if (LYZEP) {\r
261     //read the outputfile of the second run\r
262     TString outputFileName = "AnalysisResults2.root";\r
263     TString pwd(gSystem->pwd());\r
264     pwd+="/";\r
265     pwd+=outputFileName.Data();\r
266     TFile *outputFile = NULL;\r
267     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {\r
268       cout<<"WARNING: You do not have an output file:"<<endl;\r
269       cout<<"         "<<pwd.Data()<<endl;\r
270       exit(0);\r
271     } else {\r
272       outputFile = TFile::Open(pwd.Data(),"READ");\r
273     }\r
274     \r
275     // read the output file from LYZ2SUM\r
276     TString inputFileNameLYZEP = "outputLYZ2SUManalysis" ;\r
277     inputFileNameLYZEP += rptypestr;\r
278     cout<<"The input file is "<<inputFileNameLYZEP.Data()<<endl;\r
279     TFile* fInputFileLYZEP = (TFile*)outputFile->FindObjectAny(inputFileNameLYZEP.Data());\r
280     if(!fInputFileLYZEP || fInputFileLYZEP->IsZombie()) { \r
281       cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ; \r
282       break;\r
283     }\r
284     else {\r
285       TList* fInputListLYZEP = (TList*)fInputFileLYZEP->Get("LYZ2SUM");\r
286       if (!fInputListLYZEP) {cout<<"list is NULL pointer!"<<endl;}\r
287     }\r
288     cout<<"LYZEP input file/list read..."<<endl;\r
289   }\r
290   \r
291   \r
292   // Create the FMD task and add it to the manager\r
293   //===========================================================================\r
294   if (rptypestr == "FMD") {\r
295     AliFMDAnalysisTaskSE *taskfmd = NULL;\r
296     if (rptypestr == "FMD") {\r
297       taskfmd = new AliFMDAnalysisTaskSE("TaskFMD");\r
298       taskfmd->SelectCollisionCandidates(centralitysel);\r
299       mgr->AddTask(taskfmd);\r
300       \r
301       AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();\r
302       pars->Init();\r
303       pars->SetProcessPrimary(kTRUE); //for MC only\r
304       pars->SetProcessHits(kFALSE);\r
305       \r
306       //pars->SetRealData(kTRUE); //for real data\r
307       //pars->SetProcessPrimary(kFALSE); //for real data\r
308     }\r
309   }\r
310   \r
311   // Create the flow event task, add it to the manager.\r
312   //===========================================================================\r
313   AliAnalysisTaskFlowEvent *taskFE = NULL;\r
314 \r
315   if(useAfterBurner)\r
316     { \r
317       taskFE = new AliAnalysisTaskFlowEvent(Form("TaskFlowEvent %s",outputSlotName.Data()),"",doQA,1);\r
318       taskFE->SelectCollisionCandidates(centralitysel);\r
319       taskFE->SetFlow(v1,v2,v3,v4); \r
320       taskFE->SetNonFlowNumberOfTrackClones(numberOfTrackClones);\r
321       taskFE->SetAfterburnerOn();\r
322     }\r
323   else {taskFE = new AliAnalysisTaskFlowEvent(Form("TaskFlowEvent %s",outputSlotName.Data()),"",doQA); \r
324       taskFE->SelectCollisionCandidates(centralitysel);}\r
325   if (ExcludeRegion) {\r
326     taskFE->DefineDeadZone(excludeEtaMin, excludeEtaMax, excludePhiMin, excludePhiMax); \r
327   }\r
328   taskFE->SetSubeventEtaRange(minA, maxA, minB, maxB);\r
329   if (UsePhysicsSelection) {\r
330     taskFE->SelectCollisionCandidates(AliVEvent::kMB);\r
331     cout<<"Using Physics Selection"<<endl;\r
332   }\r
333   mgr->AddTask(taskFE);\r
334   \r
335   // Pass cuts for RPs and POIs to the task:\r
336   taskFE->SetCutsEvent(cutsEvent);\r
337   taskFE->SetCutsRP(cutsRP2);\r
338   taskFE->SetCutsPOI(cutsPOI);\r
339   if (cutsRP->GetParamType()==AliFlowTrackCuts::kV0)\r
340   { \r
341     //TODO: since this is set in a static object all analyses in an analysis train\r
342     //will be affected.\r
343     taskFE->SetHistWeightvsPhiMin(0.);\r
344     taskFE->SetHistWeightvsPhiMax(200.);\r
345   }\r
346 \r
347   // Create the analysis tasks, add them to the manager.\r
348   //===========================================================================\r
349   if (SP){\r
350     AliAnalysisTaskScalarProduct *taskSP = new AliAnalysisTaskScalarProduct(Form("TaskScalarProduct %s",outputSlotName.Data()),WEIGHTS[0]);\r
351     taskSP->SelectCollisionCandidates(centralitysel);\r
352     taskSP->SetRelDiffMsub(1.0);\r
353     taskSP->SetApplyCorrectionForNUA(kTRUE);\r
354     mgr->AddTask(taskSP);\r
355   }\r
356   if (LYZ1SUM){\r
357     AliAnalysisTaskLeeYangZeros *taskLYZ1SUM = new AliAnalysisTaskLeeYangZeros(Form("TaskLeeYangZerosSUM %s",outputSlotName.Data()),kTRUE);\r
358     taskLYZ1SUM->SelectCollisionCandidates(centralitysel);\r
359     taskLYZ1SUM->SetFirstRunLYZ(kTRUE);\r
360     taskLYZ1SUM->SetUseSumLYZ(kTRUE);\r
361     mgr->AddTask(taskLYZ1SUM);\r
362   }\r
363   if (LYZ1PROD){\r
364     AliAnalysisTaskLeeYangZeros *taskLYZ1PROD = new AliAnalysisTaskLeeYangZeros(Form("TaskLeeYangZerosPROD %s",outputSlotName.Data()),kTRUE);\r
365     taskLYZ1PROD->SelectCollisionCandidates(centralitysel);\r
366     taskLYZ1PROD->SetFirstRunLYZ(kTRUE);\r
367     taskLYZ1PROD->SetUseSumLYZ(kFALSE);\r
368     mgr->AddTask(taskLYZ1PROD);\r
369   }\r
370   if (LYZ2SUM){\r
371     AliAnalysisTaskLeeYangZeros *taskLYZ2SUM = new AliAnalysisTaskLeeYangZeros(Form("TaskLeeYangZerosSUM %s",outputSlotName.Data()),kFALSE);\r
372     taskLYZ2SUM->SelectCollisionCandidates(centralitysel);\r
373     taskLYZ2SUM->SetFirstRunLYZ(kFALSE);\r
374     taskLYZ2SUM->SetUseSumLYZ(kTRUE);\r
375     mgr->AddTask(taskLYZ2SUM);\r
376   }\r
377   if (LYZ2PROD){\r
378     AliAnalysisTaskLeeYangZeros *taskLYZ2PROD = new AliAnalysisTaskLeeYangZeros(Form("TaskLeeYangZerosPROD %s",outputSlotName.Data()),kFALSE);\r
379     taskLYZ2PROD->SelectCollisionCandidates(centralitysel);\r
380     taskLYZ2PROD->SetFirstRunLYZ(kFALSE);\r
381     taskLYZ2PROD->SetUseSumLYZ(kFALSE);\r
382     mgr->AddTask(taskLYZ2PROD);\r
383   }\r
384   if (LYZEP){\r
385     AliAnalysisTaskLYZEventPlane *taskLYZEP = new AliAnalysisTaskLYZEventPlane(Form("TaskLYZEventPlane %s",outputSlotName.Data()));\r
386     taskLYZEP->SelectCollisionCandidates(centralitysel);\r
387     mgr->AddTask(taskLYZEP);\r
388   }\r
389   if (GFC){\r
390     AliAnalysisTaskCumulants *taskGFC = new AliAnalysisTaskCumulants(Form("TaskCumulants %s",outputSlotName.Data()),useWeights);\r
391     taskGFC->SelectCollisionCandidates(centralitysel);\r
392     taskGFC->SetUsePhiWeights(WEIGHTS[0]); \r
393     taskGFC->SetUsePtWeights(WEIGHTS[1]);\r
394     taskGFC->SetUseEtaWeights(WEIGHTS[2]); \r
395     mgr->AddTask(taskGFC);\r
396   }\r
397   if (QC){\r
398     AliAnalysisTaskQCumulants *taskQC = new AliAnalysisTaskQCumulants(Form("TaskQCumulants %s",outputSlotName.Data()),useWeights);\r
399     taskQC->SelectCollisionCandidates(centralitysel);\r
400     taskQC->SetUsePhiWeights(WEIGHTS[0]); \r
401     taskQC->SetUsePtWeights(WEIGHTS[1]);\r
402     taskQC->SetUseEtaWeights(WEIGHTS[2]); \r
403     taskQC->SetCalculateCumulantsVsM(kFALSE);\r
404     taskQC->SetnBinsMult(10000);\r
405     taskQC->SetMinMult(0.);\r
406     taskQC->SetMaxMult(10000.);\r
407     taskQC->SetHarmonic(harmonic);\r
408     taskQC->SetApplyCorrectionForNUA(kFALSE);\r
409     taskQC->SetFillMultipleControlHistograms(kFALSE);     \r
410     mgr->AddTask(taskQC);\r
411   }\r
412   if (FQD){\r
413     AliAnalysisTaskFittingQDistribution *taskFQD = new AliAnalysisTaskFittingQDistribution(Form("TaskFittingQDistribution %s",outputSlotName.Data()),kFALSE);\r
414     taskFQD->SelectCollisionCandidates(centralitysel);\r
415     taskFQD->SetUsePhiWeights(WEIGHTS[0]); \r
416     taskFQD->SetqMin(0.);\r
417     taskFQD->SetqMax(1000.);\r
418     taskFQD->SetqNbins(10000);\r
419     mgr->AddTask(taskFQD);\r
420   }\r
421   if (MCEP){\r
422     AliAnalysisTaskMCEventPlane *taskMCEP = new AliAnalysisTaskMCEventPlane(Form("TaskMCEventPlane %s",outputSlotName.Data()));\r
423     taskMCEP->SelectCollisionCandidates(centralitysel);\r
424     mgr->AddTask(taskMCEP);\r
425   }\r
426   if (MH){\r
427     AliAnalysisTaskMixedHarmonics *taskMH = new AliAnalysisTaskMixedHarmonics(Form("TaskMixedHarmonics %s",outputSlotName.Data()),useWeights);\r
428     taskMH->SelectCollisionCandidates(centralitysel);\r
429     taskMH->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]\r
430     taskMH->SetNoOfMultipicityBins(10000);\r
431     taskMH->SetMultipicityBinWidth(1.);\r
432     taskMH->SetMinMultiplicity(1.);\r
433     taskMH->SetCorrectForDetectorEffects(kTRUE);\r
434     taskMH->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>> (Remark: two nested loops)    \r
435     taskMH->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  \r
436     mgr->AddTask(taskMH);\r
437   }  \r
438   if (NL){\r
439     AliAnalysisTaskNestedLoops *taskNL = new AliAnalysisTaskNestedLoops(Form("TaskNestedLoops %s",outputSlotName.Data()),useWeights);\r
440     taskNL->SelectCollisionCandidates(centralitysel);\r
441     taskNL->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]\r
442     taskNL->SetEvaluateNestedLoopsForRAD(kTRUE); // RAD = Relative Angle Distribution\r
443     taskNL->SetEvaluateNestedLoopsForMH(kTRUE); // evalaute <<cos[n(phi1+phi2-2phi3)]>> (Remark: three nested loops)   \r
444     taskNL->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>>  (Remark: three nested loops)   \r
445     taskNL->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  \r
446     mgr->AddTask(taskNL);\r
447   }\r
448 \r
449   // Create the output container for the data produced by the task\r
450   // Connect to the input and output containers\r
451   //===========================================================================\r
452   AliAnalysisDataContainer *cinput1 = mgr->GetCommonInputContainer();\r
453   \r
454   if (rptypestr == "FMD") {\r
455     AliAnalysisDataContainer *coutputFMD = \r
456       mgr->CreateContainer(Form("BackgroundCorrected %s",outputSlotName.Data()), TList::Class(), AliAnalysisManager::kExchangeContainer);\r
457     //input and output taskFMD     \r
458     mgr->ConnectInput(taskfmd, 0, cinput1);\r
459     mgr->ConnectOutput(taskfmd, 1, coutputFMD);\r
460     //input into taskFE\r
461     mgr->ConnectInput(taskFE,1,coutputFMD);\r
462   }\r
463   \r
464   AliAnalysisDataContainer *coutputFE = \r
465   mgr->CreateContainer(Form("FlowEventSimple %s",outputSlotName.Data()),AliFlowEventSimple::Class(),AliAnalysisManager::kExchangeContainer);\r
466   mgr->ConnectInput(taskFE,0,cinput1); \r
467   mgr->ConnectOutput(taskFE,1,coutputFE);\r
468  \r
469   if (taskFE->GetQAOn())\r
470   {\r
471     TString outputQA = fileName;\r
472     outputQA += ":QA";\r
473     AliAnalysisDataContainer* coutputFEQA = \r
474     mgr->CreateContainer(Form("QA %s",outputSlotName.Data()), TList::Class(),AliAnalysisManager::kOutputContainer,outputQA);\r
475     mgr->ConnectOutput(taskFE,2,coutputFEQA);\r
476   }\r
477 \r
478   // Create the output containers for the data produced by the analysis tasks\r
479   // Connect to the input and output containers\r
480   //===========================================================================\r
481   if (useWeights) {    \r
482     AliAnalysisDataContainer *cinputWeights = mgr->CreateContainer(Form("Weights %s",outputSlotName.Data()),\r
483                                                                    TList::Class(),AliAnalysisManager::kInputContainer); \r
484   }\r
485 \r
486   if(SP) {\r
487     TString outputSP = fileName;\r
488     outputSP += ":outputSPanalysis";\r
489     outputSP+= rptypestr;\r
490     AliAnalysisDataContainer *coutputSP = mgr->CreateContainer(Form("SP %s",outputSlotName.Data()), \r
491                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputSP); \r
492     mgr->ConnectInput(taskSP,0,coutputFE); \r
493     mgr->ConnectOutput(taskSP,1,coutputSP); \r
494     if (WEIGHTS[0]) {\r
495       mgr->ConnectInput(taskSP,1,cinputWeights);\r
496       cinputWeights->SetData(weightsList);\r
497     }\r
498   }\r
499   if(LYZ1SUM) {\r
500     TString outputLYZ1SUM = fileName;\r
501     outputLYZ1SUM += ":outputLYZ1SUManalysis";\r
502     outputLYZ1SUM+= rptypestr;\r
503     AliAnalysisDataContainer *coutputLYZ1SUM = mgr->CreateContainer(Form("LYZ1SUM %s",outputSlotName.Data()), \r
504                                                                     TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1SUM); \r
505     mgr->ConnectInput(taskLYZ1SUM,0,coutputFE);\r
506     mgr->ConnectOutput(taskLYZ1SUM,1,coutputLYZ1SUM);\r
507   }\r
508   if(LYZ1PROD) {\r
509     TString outputLYZ1PROD = fileName;\r
510     outputLYZ1PROD += ":outputLYZ1PRODanalysis";\r
511     outputLYZ1PROD+= rptypestr;\r
512     AliAnalysisDataContainer *coutputLYZ1PROD = mgr->CreateContainer(Form("LYZ1PROD %s",outputSlotName.Data()), \r
513                                                                      TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1PROD); \r
514     mgr->ConnectInput(taskLYZ1PROD,0,coutputFE); \r
515     mgr->ConnectOutput(taskLYZ1PROD,1,coutputLYZ1PROD);\r
516   }\r
517   if(LYZ2SUM) {\r
518     AliAnalysisDataContainer *cinputLYZ2SUM = mgr->CreateContainer(Form("LYZ2SUMin %s",outputSlotName.Data()),\r
519                                                                    TList::Class(),AliAnalysisManager::kInputContainer);\r
520     TString outputLYZ2SUM = fileName;\r
521     outputLYZ2SUM += ":outputLYZ2SUManalysis";\r
522     outputLYZ2SUM+= rptypestr;\r
523     \r
524     AliAnalysisDataContainer *coutputLYZ2SUM = mgr->CreateContainer(Form("LYZ2SUM %s",outputSlotName.Data()), \r
525                                                                     TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2SUM); \r
526     mgr->ConnectInput(taskLYZ2SUM,0,coutputFE); \r
527     mgr->ConnectInput(taskLYZ2SUM,1,cinputLYZ2SUM);\r
528     mgr->ConnectOutput(taskLYZ2SUM,1,coutputLYZ2SUM); \r
529     cinputLYZ2SUM->SetData(fInputListLYZ2SUM);\r
530   }\r
531   if(LYZ2PROD) {\r
532     AliAnalysisDataContainer *cinputLYZ2PROD = mgr->CreateContainer(Form("LYZ2PRODin %s",outputSlotName.Data()),\r
533                                                                     TList::Class(),AliAnalysisManager::kInputContainer);\r
534     TString outputLYZ2PROD = fileName;\r
535     outputLYZ2PROD += ":outputLYZ2PRODanalysis";\r
536     outputLYZ2PROD+= rptypestr;\r
537     \r
538     AliAnalysisDataContainer *coutputLYZ2PROD = mgr->CreateContainer(Form("LYZ2PROD %s",outputSlotName.Data()), \r
539                                                                      TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2PROD); \r
540     mgr->ConnectInput(taskLYZ2PROD,0,coutputFE); \r
541     mgr->ConnectInput(taskLYZ2PROD,1,cinputLYZ2PROD);\r
542     mgr->ConnectOutput(taskLYZ2PROD,1,coutputLYZ2PROD); \r
543     cinputLYZ2PROD->SetData(fInputListLYZ2PROD);\r
544   }\r
545   if(LYZEP) {\r
546     AliAnalysisDataContainer *cinputLYZEP = mgr->CreateContainer(Form("LYZEPin %s",outputSlotName.Data()),\r
547                                                                  TList::Class(),AliAnalysisManager::kInputContainer);\r
548     TString outputLYZEP = fileName;\r
549     outputLYZEP += ":outputLYZEPanalysis";\r
550     outputLYZEP+= rptypestr;\r
551     \r
552     AliAnalysisDataContainer *coutputLYZEP = mgr->CreateContainer(Form("LYZEP %s",outputSlotName.Data()), \r
553                                                                   TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZEP); \r
554     mgr->ConnectInput(taskLYZEP,0,coutputFE); \r
555     mgr->ConnectInput(taskLYZEP,1,cinputLYZEP);\r
556     mgr->ConnectOutput(taskLYZEP,1,coutputLYZEP); \r
557     cinputLYZEP->SetData(fInputListLYZEP);\r
558   }\r
559   if(GFC) {\r
560     TString outputGFC = fileName;\r
561     outputGFC += ":outputGFCanalysis";\r
562     outputGFC+= rptypestr;\r
563     \r
564     AliAnalysisDataContainer *coutputGFC = mgr->CreateContainer(Form("GFC %s",outputSlotName.Data()), \r
565                                                                 TList::Class(),AliAnalysisManager::kOutputContainer,outputGFC); \r
566     mgr->ConnectInput(taskGFC,0,coutputFE); \r
567     mgr->ConnectOutput(taskGFC,1,coutputGFC);\r
568     if (useWeights) {\r
569       mgr->ConnectInput(taskGFC,1,cinputWeights);\r
570       cinputWeights->SetData(weightsList);\r
571     } \r
572   }\r
573   if(QC) {\r
574     TString outputQC = fileName;\r
575     outputQC += ":outputQCanalysis";\r
576     outputQC+= rptypestr;\r
577 \r
578     AliAnalysisDataContainer *coutputQC = mgr->CreateContainer(Form("QC %s",outputSlotName.Data()), \r
579                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputQC); \r
580     mgr->ConnectInput(taskQC,0,coutputFE); \r
581     mgr->ConnectOutput(taskQC,1,coutputQC);\r
582     if (useWeights) {\r
583       mgr->ConnectInput(taskQC,1,cinputWeights);\r
584       cinputWeights->SetData(weightsList);\r
585     }\r
586   }\r
587   if(FQD) {\r
588     TString outputFQD = fileName;\r
589     outputFQD += ":outputFQDanalysis";\r
590     outputFQD+= rptypestr;\r
591     \r
592     AliAnalysisDataContainer *coutputFQD = mgr->CreateContainer(Form("FQD %s",outputSlotName.Data()), \r
593                                                                 TList::Class(),AliAnalysisManager::kOutputContainer,outputFQD); \r
594     mgr->ConnectInput(taskFQD,0,coutputFE); \r
595     mgr->ConnectOutput(taskFQD,1,coutputFQD);\r
596     if(useWeights) {\r
597       mgr->ConnectInput(taskFQD,1,cinputWeights);\r
598       cinputWeights->SetData(weightsList);\r
599     } \r
600   }\r
601   if(MCEP) {\r
602     TString outputMCEP = fileName;\r
603     outputMCEP += ":outputMCEPanalysis";\r
604     outputMCEP+= rptypestr;\r
605     \r
606     AliAnalysisDataContainer *coutputMCEP = mgr->CreateContainer(Form("MCEP %s",outputSlotName.Data()), \r
607                                                                  TList::Class(),AliAnalysisManager::kOutputContainer,outputMCEP); \r
608     mgr->ConnectInput(taskMCEP,0,coutputFE);\r
609     mgr->ConnectOutput(taskMCEP,1,coutputMCEP); \r
610   }\r
611   if(MH) {\r
612     TString outputMH = fileName;\r
613     outputMH += ":outputMHanalysis";\r
614     outputMH += rptypestr;\r
615         \r
616     AliAnalysisDataContainer *coutputMH = mgr->CreateContainer(Form("MH %s",outputSlotName.Data()), \r
617                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputMH); \r
618     mgr->ConnectInput(taskMH,0,coutputFE); \r
619     mgr->ConnectOutput(taskMH,1,coutputMH); \r
620     //if (useWeights) {\r
621     //  mgr->ConnectInput(taskMH,1,cinputWeights);\r
622     //  cinputWeights->SetData(weightsList);\r
623     //} \r
624   }\r
625   if(NL) {\r
626     TString outputNL = fileName;\r
627     outputNL += ":outputNLanalysis";\r
628     outputNL += rptypestr;\r
629 \r
630     AliAnalysisDataContainer *coutputNL = mgr->CreateContainer(Form("NL %s",outputSlotName.Data()), \r
631                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputNL); \r
632     mgr->ConnectInput(taskNL,0,coutputFE);\r
633     mgr->ConnectOutput(taskNL,1,coutputNL);\r
634     //if (useWeights) {\r
635     //  mgr->ConnectInput(taskNL,1,cinputWeights);\r
636     //  cinputWeights->SetData(weightsList);\r
637     //} \r
638   }\r
639 \r
640   ///////////////////////////////////////////////////////////////////////////////////////////\r
641   if (runQAtask)\r
642   {\r
643     AliAnalysisTaskQAflow* taskQAflow = new AliAnalysisTaskQAflow(Form("TaskQAflow %s",outputSlotName.Data()));\r
644     taskQAflow->SelectCollisionCandidates(centralitysel);\r
645     taskQAflow->SetEventCuts(cutsEvent);\r
646     taskQAflow->SetTrackCuts(cutsRP);\r
647     taskQAflow->SetFillNTuple(FillQAntuple);\r
648     taskQAflow->SetDoCorrelations(DoQAcorrelations);\r
649     mgr->AddTask(taskQAflow);\r
650     \r
651     Printf("outputSlotName %s",outputSlotName.Data());\r
652     TString taskQAoutputFileName(fileNameBase);\r
653     taskQAoutputFileName.Append("_QA.root");\r
654     AliAnalysisDataContainer* coutputQAtask = mgr->CreateContainer(Form("flowQA %s",outputSlotName.Data()),\r
655                                               TObjArray::Class(),\r
656                                               AliAnalysisManager::kOutputContainer,\r
657                                               taskQAoutputFileName);\r
658     AliAnalysisDataContainer* coutputQAtaskTree = mgr->CreateContainer(Form("flowQAntuple %s",outputSlotName.Data()),\r
659                                               TNtuple::Class(),\r
660                                               AliAnalysisManager::kOutputContainer,\r
661                                               taskQAoutputFileName);\r
662     mgr->ConnectInput(taskQAflow,0,mgr->GetCommonInputContainer());\r
663     mgr->ConnectInput(taskQAflow,1,coutputFE);\r
664     mgr->ConnectOutput(taskQAflow,1,coutputQAtask);\r
665     if (FillQAntuple) mgr->ConnectOutput(taskQAflow,2,coutputQAtaskTree);\r
666   }\r
667 }\r
668 \r
669 \r
670 \r
671 \r
672 \r