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