]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/macros/AddTaskFlow.C
update for LYZ
[u/mrichter/AliRoot.git] / PWG2 / FLOW / macros / AddTaskFlow.C
1 /////////////////////////////////////////////////////////////////////////////////////////////
2 //
3 // AddTask* macro for flow analysis
4 // Creates a Flow Event task and adds it to the analysis manager.
5 // Sets the cuts using the correction framework (CORRFW) classes.
6 // Also creates Flow Analysis tasks and connects them to the output of the flow event task.
7 //
8 /////////////////////////////////////////////////////////////////////////////////////////////
9
10
11 // SETTING THE CUTS
12
13 // event selection
14 const Int_t multminESD = 3;  //used for CORRFW cuts 
15 //const Int_t multmaxESD = 1000; //used for CORRFW cuts 
16 const Int_t multmaxESD = 1000000; //used for CORRFW cuts 
17
18 const Int_t multmin = 3;     //used for AliFlowEventSimple (to set the centrality)
19 //const Int_t multmax = 100;     //used for AliFlowEventSimple (to set the centrality)
20 const Int_t multmax = 1000000;     //used for AliFlowEventSimple (to set the centrality)
21
22
23 // For RP selection
24 const Double_t ptmin1 = 0.0;
25 const Double_t ptmax1 = 10.0;
26 const Double_t ymin1  = -1.;
27 const Double_t ymax1  = 1.;
28 const Int_t mintrackrefsTPC1 = 2; //2;
29 const Int_t mintrackrefsITS1 = 3;
30 const Int_t charge1 = 1;
31 Bool_t UsePIDforRP = kFALSE;
32 const Int_t PDG1 = 211;
33 const Int_t minclustersTPC1 = 70;//50;
34 const Int_t maxnsigmatovertex1 = 3;
35
36 // For for POI selection
37 const Double_t ptmin2 = 0.0;
38 const Double_t ptmax2 = 10.0;
39 const Double_t ymin2  = -1.;
40 const Double_t ymax2  = 1.;
41 const Int_t mintrackrefsTPC2 = 2; //2;
42 const Int_t mintrackrefsITS2 = 3;
43 const Int_t charge2 = 1;
44 Bool_t UsePIDforPOI = kFALSE;
45 const Int_t PDG2 = 321;
46 const Int_t minclustersTPC2 = 70; //50;
47 const Int_t maxnsigmatovertex2 = 3;
48
49 // For manipulating the event (for testing purposes)
50 const Bool_t AddToEvent = kFALSE;
51 Double_t ellipticflow = 0.05;
52
53
54 AliAnalysisTaskFlowEvent* AddTaskFlow(TString type, Bool_t* METHODS, Bool_t QA, Bool_t* WEIGHTS)
55 {
56   //boleans for the methods
57   Bool_t SP       = METHODS[0];
58   Bool_t LYZ1SUM  = METHODS[1];
59   Bool_t LYZ1PROD = METHODS[2];
60   Bool_t LYZ2SUM  = METHODS[3];
61   Bool_t LYZ2PROD = METHODS[4];
62   Bool_t LYZEP    = METHODS[5];
63   Bool_t GFC      = METHODS[6];
64   Bool_t QC       = METHODS[7];
65   Bool_t FQD      = METHODS[8];
66   Bool_t MCEP     = METHODS[9];   
67  
68   //for using weights
69   Bool_t useWeights  = WEIGHTS[0] || WEIGHTS[1] || WEIGHTS[2];
70   if (useWeights) cout<<"Weights are used"<<endl;
71   else cout<<"Weights are not used"<<endl;
72
73
74   // Get the pointer to the existing analysis manager via the static access method.
75   //==============================================================================
76   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
77   if (!mgr) {
78     Error("AddTaskFlowEvent", "No analysis manager to connect to.");
79     return NULL;
80   }   
81   
82   // Check the analysis type using the event handlers connected to the analysis
83   // manager. The availability of MC handler cann also be checked here.
84   //==============================================================================
85   if (!mgr->GetInputEventHandler()) {
86     ::Error("AddTaskFlowEvent", "This task requires an input event handler");
87     return NULL;
88   }  
89     
90   // Open external input files
91   //===========================================================================
92   //weights: 
93   TFile *weightsFile = NULL;
94   TList *weightsList = NULL;
95
96   if(useWeights) {
97     //open the file with the weights:
98     weightsFile = TFile::Open("weights.root","READ");
99     if(weightsFile) {
100       //access the list which holds the histos with weigths:
101       weightsList = (TList*)weightsFile->Get("weights");
102     }
103     else {
104       cout<<" WARNING: the file <weights.root> with weights from the previous run was not available."<<endl;
105       break;
106     } 
107   }
108     
109   if (LYZ2SUM || LYZ2PROD) {
110     //read the outputfile of the first run
111     TString outputFileName = "AnalysisResults1.root"; 
112     TString pwd(gSystem->pwd());
113     pwd+="/";
114     pwd+=outputFileName.Data();
115     TFile *outputFile = NULL;
116     if(gSystem->AccessPathName(pwd.Data(),kFileExists))
117       {
118         cout<<"WARNING: You do not have an output file:"<<endl;
119         cout<<"         "<<pwd.Data()<<endl;
120         exit(0);
121       } else 
122         {
123           outputFile = TFile::Open(pwd.Data(),"READ");
124         }
125  
126     if (LYZ2SUM){  
127       // read the output file from LYZ1SUM 
128       TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis" ;
129       inputFileNameLYZ2SUM += type;
130       //inputFileNameLYZ2SUM += ".root";
131       cout<<"The input file is "<<inputFileNameLYZ2SUM.Data()<<endl;
132       TFile* fInputFileLYZ2SUM = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2SUM.Data());
133       
134       if(!fInputFileLYZ2SUM || fInputFileLYZ2SUM->IsZombie()) { 
135         cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ; 
136         break;
137       }
138       else {
139         TList* fInputListLYZ2SUM = (TList*)fInputFileLYZ2SUM->Get("cobjLYZ1SUM");
140         if (!fInputListLYZ2SUM) {cout<<"list is NULL pointer!"<<endl;}
141       }
142       cout<<"LYZ2SUM input file/list read..."<<endl;
143     }
144
145     if (LYZ2PROD){  
146       // read the output file from LYZ1PROD 
147       TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis" ;
148       inputFileNameLYZ2PROD += type;
149       //inputFileNameLYZ2PROD += ".root";
150       cout<<"The input file is "<<inputFileNameLYZ2PROD.Data()<<endl;
151       TFile* fInputFileLYZ2PROD =  (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2PROD.Data());
152       if(!fInputFileLYZ2PROD || fInputFileLYZ2PROD->IsZombie()) { 
153         cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ; 
154         break;
155       }
156       else {
157         TList* fInputListLYZ2PROD = (TList*)fInputFileLYZ2PROD->Get("cobjLYZ1PROD");
158         if (!fInputListLYZ2PROD) {cout<<"list is NULL pointer!"<<endl;}
159       }
160       cout<<"LYZ2PROD input file/list read..."<<endl;
161     }
162   }
163   
164
165   if (LYZEP) {
166     // read the output file from LYZ2SUM
167     TString outputFileName = "AnalysisResults2.root"; 
168     TString pwd(gSystem->pwd());
169     pwd+="/";
170     pwd+=outputFileName.Data();
171     TFile *outputFile = NULL;
172     if(gSystem->AccessPathName(pwd.Data(),kFileExists))
173       {
174         cout<<"WARNING: You do not have an output file:"<<endl;
175         cout<<"         "<<pwd.Data()<<endl;
176         exit(0);
177       } else 
178         {
179           outputFile = TFile::Open(pwd.Data(),"READ");
180         }
181
182     TString inputFileNameLYZEP = "outputLYZ2SUManalysis" ;
183     inputFileNameLYZEP += type;
184     //inputFileNameLYZEP += ".root";
185     cout<<"The input file is "<<inputFileNameLYZEP.Data()<<endl;
186     TFile* fInputFileLYZEP = (TFile*)outputFile->FindObjectAny(inputFileNameLYZEP.Data());
187     if(!fInputFileLYZEP || fInputFileLYZEP->IsZombie()) { 
188       cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ; 
189       break;
190     }
191     else {
192       TList* fInputListLYZEP = (TList*)fInputFileLYZEP->Get("cobjLYZ2SUM");
193       if (!fInputListLYZEP) {cout<<"list is NULL pointer!"<<endl;}
194     }
195     cout<<"LYZEP input file/list read..."<<endl;
196   }
197   
198
199
200   // Create the task, add it to the manager.
201   //===========================================================================
202   AliAnalysisTaskFlowEvent *taskFE = NULL;
203   if (QA) { 
204     if(AddToEvent) { 
205       taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",kTRUE,1);
206       taskFE->SetEllipticFlowValue(ellipticflow); }    //TEST
207     else {taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",kTRUE); }
208     taskFE->SetAnalysisType(type);
209     taskFE->SetMinMult(multmin);
210     taskFE->SetMaxMult(multmax);
211     mgr->AddTask(taskFE);
212   }
213   else { 
214     taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",kFALSE); 
215     taskFE->SetAnalysisType(type);
216     taskFE->SetMinMult(multmin);
217     taskFE->SetMaxMult(multmax);
218     mgr->AddTask(taskFE);
219   }
220  
221   // Create cuts using the correction framework (CORRFW)
222   //===========================================================================
223   if (QA){
224     //Set TList for the QA histograms
225     TList* qaRP  = new TList(); 
226     TList* qaPOI = new TList();
227   }
228
229   //############# event cuts on multiplicity
230   AliCFEventGenCuts* mcEventCuts = new AliCFEventGenCuts("mcEventCuts","MC-level event cuts");
231   mcEventCuts->SetNTracksCut(multminESD,multmaxESD); 
232   if (QA) { 
233     mcEventCuts->SetQAOn(qaRP);
234   }
235   AliCFEventRecCuts* recEventCuts = new AliCFEventRecCuts("recEventCuts","rec-level event cuts");
236   recEventCuts->SetNTracksCut(multminESD,multmaxESD); 
237   if (QA) { 
238     recEventCuts->SetQAOn(qaRP);
239   }
240
241   //############# cuts on MC
242   AliCFTrackKineCuts* mcKineCuts1 = new AliCFTrackKineCuts("mcKineCuts1","MC-level kinematic cuts");
243   mcKineCuts1->SetPtRange(ptmin1,ptmax1);
244   mcKineCuts1->SetRapidityRange(ymin1,ymax1);
245   //mcKineCuts1->SetChargeMC(charge1);
246   if (QA) { 
247     mcKineCuts1->SetQAOn(qaRP);
248   }
249
250   AliCFTrackKineCuts* mcKineCuts2 = new AliCFTrackKineCuts("mcKineCuts2","MC-level kinematic cuts");
251   mcKineCuts2->SetPtRange(ptmin2,ptmax2);
252   mcKineCuts2->SetRapidityRange(ymin2,ymax2);
253   //mcKineCuts2->SetChargeMC(charge2);
254   if (QA) { 
255     mcKineCuts2->SetQAOn(qaPOI);
256   }
257   
258   AliCFParticleGenCuts* mcGenCuts1 = new AliCFParticleGenCuts("mcGenCuts1","MC particle generation cuts for RP");
259   mcGenCuts1->SetRequireIsPrimary();
260   if (UsePIDforRP) {mcGenCuts1->SetRequirePdgCode(PDG1);}
261   if (QA) { 
262     mcGenCuts1->SetQAOn(qaRP);
263   }
264   
265   AliCFParticleGenCuts* mcGenCuts2 = new AliCFParticleGenCuts("mcGenCuts2","MC particle generation cuts for POI");
266   mcGenCuts2->SetRequireIsPrimary();
267   if (UsePIDforPOI) {mcGenCuts2->SetRequirePdgCode(PDG2);}
268   if (QA) { 
269     mcGenCuts2->SetQAOn(qaPOI);
270   }
271
272   //############# Acceptance Cuts  
273   AliCFAcceptanceCuts *mcAccCuts1 = new AliCFAcceptanceCuts("mcAccCuts1","MC acceptance cuts");
274   mcAccCuts1->SetMinNHitITS(mintrackrefsITS1);
275   mcAccCuts1->SetMinNHitTPC(mintrackrefsTPC1);
276   if (QA) { 
277     mcAccCuts1->SetQAOn(qaRP);
278   }
279   
280   AliCFAcceptanceCuts *mcAccCuts2 = new AliCFAcceptanceCuts("mcAccCuts2","MC acceptance cuts");
281   mcAccCuts2->SetMinNHitITS(mintrackrefsITS2);
282   mcAccCuts2->SetMinNHitTPC(mintrackrefsTPC2);
283   if (QA) { 
284     mcAccCuts2->SetQAOn(qaPOI);
285   }
286   //############# Rec-Level kinematic cuts
287   AliCFTrackKineCuts *recKineCuts1 = new AliCFTrackKineCuts("recKineCuts1","rec-level kine cuts");
288   recKineCuts1->SetPtRange(ptmin1,ptmax1);
289   recKineCuts1->SetRapidityRange(ymin1,ymax1);
290   //recKineCuts1->SetChargeRec(charge1);
291   if (QA) { 
292     recKineCuts1->SetQAOn(qaRP);
293   }
294   
295   AliCFTrackKineCuts *recKineCuts2 = new AliCFTrackKineCuts("recKineCuts2","rec-level kine cuts");
296   recKineCuts2->SetPtRange(ptmin2,ptmax2);
297   recKineCuts2->SetRapidityRange(ymin2,ymax2);
298   //recKineCuts2->SetChargeRec(charge2);
299   if (QA) { 
300     recKineCuts2->SetQAOn(qaPOI);
301   }
302   
303   AliCFTrackQualityCuts *recQualityCuts1 = new AliCFTrackQualityCuts("recQualityCuts1","rec-level quality cuts");
304   recQualityCuts1->SetMinNClusterTPC(minclustersTPC1);
305   recQualityCuts1->SetStatus(AliESDtrack::kITSrefit);
306   if (QA) { 
307     recQualityCuts1->SetQAOn(qaRP);
308   }
309   AliCFTrackQualityCuts *recQualityCuts2 = new AliCFTrackQualityCuts("recQualityCuts2","rec-level quality cuts");
310   recQualityCuts2->SetMinNClusterTPC(minclustersTPC2);
311   recQualityCuts2->SetStatus(AliESDtrack::kITSrefit);
312   if (QA) { 
313     recQualityCuts2->SetQAOn(qaPOI);
314   }
315   
316   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts1 = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts1","rec-level isPrimary cuts");
317   recIsPrimaryCuts1->SetMaxNSigmaToVertex(maxnsigmatovertex1);
318   if (QA) { 
319     recIsPrimaryCuts1->SetQAOn(qaRP);
320   }
321   
322   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts2 = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts2","rec-level isPrimary cuts");
323   recIsPrimaryCuts2->SetMaxNSigmaToVertex(maxnsigmatovertex2);
324   if (QA) { 
325     recIsPrimaryCuts2->SetQAOn(qaPOI);
326   }
327   
328   int n_species = AliPID::kSPECIES ;
329   Double_t* prior = new Double_t[n_species];
330   
331   prior[0] = 0.0244519 ;
332   prior[1] = 0.0143988 ;
333   prior[2] = 0.805747  ;
334   prior[3] = 0.0928785 ;
335   prior[4] = 0.0625243 ;
336   
337   AliCFTrackCutPid* cutPID1 = NULL;
338   if(UsePIDforRP) {
339     AliCFTrackCutPid* cutPID1 = new AliCFTrackCutPid("cutPID1","ESD_PID for RP") ;
340     cutPID1->SetPriors(prior);
341     cutPID1->SetProbabilityCut(0.0);
342     cutPID1->SetDetectors("TPC TOF");
343     switch(TMath::Abs(PDG1)) {
344     case 11   : cutPID1->SetParticleType(AliPID::kElectron, kTRUE); break;
345     case 13   : cutPID1->SetParticleType(AliPID::kMuon    , kTRUE); break;
346     case 211  : cutPID1->SetParticleType(AliPID::kPion    , kTRUE); break;
347     case 321  : cutPID1->SetParticleType(AliPID::kKaon    , kTRUE); break;
348     case 2212 : cutPID1->SetParticleType(AliPID::kProton  , kTRUE); break;
349     default   : printf("UNDEFINED PID\n"); break;
350     }
351     if (QA) { 
352       cutPID1->SetQAOn(qaRP); 
353     }
354   }
355   
356   AliCFTrackCutPid* cutPID2 = NULL;
357   if (UsePIDforPOI) {
358     AliCFTrackCutPid* cutPID2 = new AliCFTrackCutPid("cutPID2","ESD_PID for POI") ;
359     cutPID2->SetPriors(prior);
360     cutPID2->SetProbabilityCut(0.0);
361     cutPID2->SetDetectors("TPC TOF");
362     switch(TMath::Abs(PDG2)) {
363     case 11   : cutPID2->SetParticleType(AliPID::kElectron, kTRUE); break;
364     case 13   : cutPID2->SetParticleType(AliPID::kMuon    , kTRUE); break;
365     case 211  : cutPID2->SetParticleType(AliPID::kPion    , kTRUE); break;
366     case 321  : cutPID2->SetParticleType(AliPID::kKaon    , kTRUE); break;
367     case 2212 : cutPID2->SetParticleType(AliPID::kProton  , kTRUE); break;
368     default   : printf("UNDEFINED PID\n"); break;
369     }
370     if (QA) { 
371       cutPID2->SetQAOn(qaPOI);
372     }
373   }
374   
375   printf("CREATE EVENT CUTS\n");
376   TObjArray* mcEventList = new TObjArray(0);
377   mcEventList->AddLast(mcEventCuts);
378     
379   TObjArray* recEventList = new TObjArray(0);
380   recEventList->AddLast(recEventCuts);
381
382   printf("CREATE MC KINE CUTS\n");
383   TObjArray* mcList1 = new TObjArray(0);
384   mcList1->AddLast(mcKineCuts1);
385   mcList1->AddLast(mcGenCuts1);
386   
387   TObjArray* mcList2 = new TObjArray(0);
388   mcList2->AddLast(mcKineCuts2);
389   mcList2->AddLast(mcGenCuts2);
390   
391   printf("CREATE ACCEPTANCE CUTS\n");
392   TObjArray* accList1 = new TObjArray(0) ;
393   accList1->AddLast(mcAccCuts1);
394   
395   TObjArray* accList2 = new TObjArray(0) ;
396   accList2->AddLast(mcAccCuts2);
397   
398   printf("CREATE RECONSTRUCTION CUTS\n");
399   TObjArray* recList1 = new TObjArray(0) ;
400   recList1->AddLast(recKineCuts1);
401   recList1->AddLast(recQualityCuts1);
402   recList1->AddLast(recIsPrimaryCuts1);
403   
404   TObjArray* recList2 = new TObjArray(0) ;
405   recList2->AddLast(recKineCuts2);
406   recList2->AddLast(recQualityCuts2);
407   recList2->AddLast(recIsPrimaryCuts2);
408   
409   printf("CREATE PID CUTS\n");
410   TObjArray* fPIDCutList1 = new TObjArray(0) ;
411   if(UsePIDforRP) {fPIDCutList1->AddLast(cutPID1);}
412   
413   TObjArray* fPIDCutList2 = new TObjArray(0) ;
414   if (UsePIDforPOI)  {fPIDCutList2->AddLast(cutPID2);}
415   
416   printf("CREATE INTERFACE AND CUTS\n");
417   AliCFManager* cfmgr1 = new AliCFManager();
418   cfmgr1->SetNStepEvent(3);
419   cfmgr1->SetEventCutsList(AliCFManager::kEvtGenCuts,mcEventList); 
420   cfmgr1->SetEventCutsList(AliCFManager::kEvtRecCuts,recEventList); 
421   cfmgr1->SetNStepParticle(4); 
422   cfmgr1->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList1);
423   cfmgr1->SetParticleCutsList(AliCFManager::kPartAccCuts,accList1);
424   cfmgr1->SetParticleCutsList(AliCFManager::kPartRecCuts,recList1);
425   cfmgr1->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList1);
426   
427   AliCFManager* cfmgr2 = new AliCFManager();
428   cfmgr1->SetNStepEvent(3);
429   cfmgr1->SetEventCutsList(AliCFManager::kEvtGenCuts,mcEventList); 
430   cfmgr1->SetEventCutsList(AliCFManager::kEvtRecCuts,recEventList); 
431   cfmgr2->SetNStepParticle(4); 
432   cfmgr2->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList2);
433   cfmgr2->SetParticleCutsList(AliCFManager::kPartAccCuts,accList2);
434   cfmgr2->SetParticleCutsList(AliCFManager::kPartRecCuts,recList2);
435   cfmgr2->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList2);
436   
437   if (QA) {
438     taskFE->SetQAList1(qaRP);
439     taskFE->SetQAList2(qaPOI);
440   }
441   taskFE->SetCFManager1(cfmgr1);
442   taskFE->SetCFManager2(cfmgr2);
443
444
445
446   // Create the analysis tasks, add them to the manager.
447   //===========================================================================
448   if (SP){
449     AliAnalysisTaskScalarProduct *taskSP = new AliAnalysisTaskScalarProduct("TaskScalarProduct");
450     mgr->AddTask(taskSP);
451   }
452   if (LYZ1SUM){
453     AliAnalysisTaskLeeYangZeros *taskLYZ1SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kTRUE);
454     taskLYZ1SUM->SetFirstRunLYZ(kTRUE);
455     taskLYZ1SUM->SetUseSumLYZ(kTRUE);
456     mgr->AddTask(taskLYZ1SUM);
457   }
458   if (LYZ1PROD){
459     AliAnalysisTaskLeeYangZeros *taskLYZ1PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kTRUE);
460     taskLYZ1PROD->SetFirstRunLYZ(kTRUE);
461     taskLYZ1PROD->SetUseSumLYZ(kFALSE);
462     mgr->AddTask(taskLYZ1PROD);
463   }
464   if (LYZ2SUM){
465     AliAnalysisTaskLeeYangZeros *taskLYZ2SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kFALSE);
466     taskLYZ2SUM->SetFirstRunLYZ(kFALSE);
467     taskLYZ2SUM->SetUseSumLYZ(kTRUE);
468     mgr->AddTask(taskLYZ2SUM);
469   }
470   if (LYZ2PROD){
471     AliAnalysisTaskLeeYangZeros *taskLYZ2PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kFALSE);
472     taskLYZ2PROD->SetFirstRunLYZ(kFALSE);
473     taskLYZ2PROD->SetUseSumLYZ(kFALSE);
474     mgr->AddTask(taskLYZ2PROD);
475   }
476   if (LYZEP){
477     AliAnalysisTaskLYZEventPlane *taskLYZEP = new AliAnalysisTaskLYZEventPlane("TaskLYZEventPlane");
478     mgr->AddTask(taskLYZEP);
479   }
480   if (GFC){
481     AliAnalysisTaskCumulants *taskGFC = new AliAnalysisTaskCumulants("TaskCumulants",useWeights);
482     taskGFC->SetUsePhiWeights(WEIGHTS[0]); 
483     taskGFC->SetUsePtWeights(WEIGHTS[1]);
484     taskGFC->SetUseEtaWeights(WEIGHTS[2]); 
485     mgr->AddTask(taskGFC);
486   }
487   if (QC){
488     AliAnalysisTaskQCumulants *taskQC = new AliAnalysisTaskQCumulants("TaskQCumulants",useWeights);
489     taskQC->SetUsePhiWeights(WEIGHTS[0]); 
490     taskQC->SetUsePtWeights(WEIGHTS[1]);
491     taskQC->SetUseEtaWeights(WEIGHTS[2]); 
492     mgr->AddTask(taskQC);
493   }
494   if (FQD){
495     AliAnalysisTaskFittingQDistribution *taskFQD = new AliAnalysisTaskFittingQDistribution("TaskFittingQDistribution",kFALSE);
496     taskFQD->SetUsePhiWeights(WEIGHTS[0]); 
497     mgr->AddTask(taskFQD);
498   }
499   if (MCEP){
500     AliAnalysisTaskMCEventPlane *taskMCEP = new AliAnalysisTaskMCEventPlane("TaskMCEventPlane");
501     mgr->AddTask(taskMCEP);
502   }
503
504   // Create the output container for the data produced by the task
505   // Connect to the input and output containers
506   //===========================================================================
507   AliAnalysisDataContainer *cinput1 = mgr->GetCommonInputContainer();
508   AliAnalysisDataContainer *coutputFE = mgr->CreateContainer("cobjFlowEventSimple",  AliFlowEventSimple::Class(),AliAnalysisManager::kExchangeContainer);
509   mgr->ConnectInput(taskFE,0,cinput1); 
510   mgr->ConnectOutput(taskFE,0,coutputFE);
511
512   if (QA) { 
513     TString qaNameRPFE = AliAnalysisManager::GetCommonFileName();
514     qaNameRPFE += ":QAforRP_FE_";
515     qaNameRPFE += type;
516
517     AliAnalysisDataContainer *coutputQA1FE =
518       mgr->CreateContainer("QARPFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNameRPFE); 
519     
520     TString qaNamePOIFE = AliAnalysisManager::GetCommonFileName();
521     qaNamePOIFE += ":QAforPOI_FE_";
522     qaNamePOIFE += type;
523         
524     AliAnalysisDataContainer *coutputQA2FE =
525       mgr->CreateContainer("QAPOIFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNamePOIFE); 
526
527     mgr->ConnectOutput(taskFE,1,coutputQA1FE); 
528     mgr->ConnectOutput(taskFE,2,coutputQA2FE); 
529   }
530
531   // Create the output containers for the data produced by the analysis tasks
532   // Connect to the input and output containers
533   //===========================================================================
534   if (useWeights) {    
535     AliAnalysisDataContainer *cinputWeights = mgr->CreateContainer("cobjWeights",TList::Class(),AliAnalysisManager::kInputContainer); 
536   }
537
538   if(SP) {
539     TString outputSP = AliAnalysisManager::GetCommonFileName();
540     outputSP += ":outputSPanalysis";
541     outputSP+= type;
542     
543     AliAnalysisDataContainer *coutputSP = mgr->CreateContainer("cobjSP", TList::Class(),AliAnalysisManager::kOutputContainer,outputSP); 
544     mgr->ConnectInput(taskSP,0,coutputFE); 
545     mgr->ConnectOutput(taskSP,0,coutputSP); 
546   }
547   if(LYZ1SUM) {
548     TString outputLYZ1SUM = AliAnalysisManager::GetCommonFileName();
549     outputLYZ1SUM += ":outputLYZ1SUManalysis";
550     outputLYZ1SUM+= type;
551     
552     AliAnalysisDataContainer *coutputLYZ1SUM = mgr->CreateContainer("cobjLYZ1SUM", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1SUM); 
553     mgr->ConnectInput(taskLYZ1SUM,0,coutputFE); 
554     mgr->ConnectOutput(taskLYZ1SUM,0,coutputLYZ1SUM); 
555   }
556   if(LYZ1PROD) {
557     TString outputLYZ1PROD = AliAnalysisManager::GetCommonFileName();
558     outputLYZ1PROD += ":outputLYZ1PRODanalysis";
559     outputLYZ1PROD+= type;
560     
561     AliAnalysisDataContainer *coutputLYZ1PROD = mgr->CreateContainer("cobjLYZ1PROD", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1PROD); 
562     mgr->ConnectInput(taskLYZ1PROD,0,coutputFE); 
563     mgr->ConnectOutput(taskLYZ1PROD,0,coutputLYZ1PROD);
564   }
565   if(LYZ2SUM) {
566     AliAnalysisDataContainer *cinputLYZ2SUM = mgr->CreateContainer("cobjLYZ2SUMin",TList::Class(),AliAnalysisManager::kInputContainer);
567     TString outputLYZ2SUM = AliAnalysisManager::GetCommonFileName();
568     outputLYZ2SUM += ":outputLYZ2SUManalysis";
569     outputLYZ2SUM+= type;
570     
571     AliAnalysisDataContainer *coutputLYZ2SUM = mgr->CreateContainer("cobjLYZ2SUM", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2SUM); 
572     mgr->ConnectInput(taskLYZ2SUM,0,coutputFE); 
573     mgr->ConnectInput(taskLYZ2SUM,1,cinputLYZ2SUM);
574     mgr->ConnectOutput(taskLYZ2SUM,0,coutputLYZ2SUM); 
575     cinputLYZ2SUM->SetData(fInputListLYZ2SUM);
576   }
577   if(LYZ2PROD) {
578     AliAnalysisDataContainer *cinputLYZ2PROD = mgr->CreateContainer("cobjLYZ2PRODin",TList::Class(),AliAnalysisManager::kInputContainer);
579     TString outputLYZ2PROD = AliAnalysisManager::GetCommonFileName();
580     outputLYZ2PROD += ":outputLYZ2PRODanalysis";
581     outputLYZ2PROD+= type;
582     
583     AliAnalysisDataContainer *coutputLYZ2PROD = mgr->CreateContainer("cobjLYZ2PROD", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2PROD); 
584     mgr->ConnectInput(taskLYZ2PROD,0,coutputFE); 
585     mgr->ConnectInput(taskLYZ2PROD,1,cinputLYZ2PROD);
586     mgr->ConnectOutput(taskLYZ2PROD,0,coutputLYZ2PROD); 
587     cinputLYZ2PROD->SetData(fInputListLYZ2PROD);
588   }
589   if(LYZEP) {
590     AliAnalysisDataContainer *cinputLYZEP = mgr->CreateContainer("cobjLYZEPin",TList::Class(),AliAnalysisManager::kInputContainer);
591     TString outputLYZEP = AliAnalysisManager::GetCommonFileName();
592     outputLYZEP += ":outputLYZEPanalysis";
593     outputLYZEP+= type;
594     
595     AliAnalysisDataContainer *coutputLYZEP = mgr->CreateContainer("cobjLYZEP", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZEP); 
596     mgr->ConnectInput(taskLYZEP,0,coutputFE); 
597     mgr->ConnectInput(taskLYZEP,1,cinputLYZEP);
598     mgr->ConnectOutput(taskLYZEP,0,coutputLYZEP); 
599     cinputLYZEP->SetData(fInputListLYZEP);
600   }
601   if(GFC) {
602     TString outputGFC = AliAnalysisManager::GetCommonFileName();
603     outputGFC += ":outputGFCanalysis";
604     outputGFC+= type;
605     
606     AliAnalysisDataContainer *coutputGFC = mgr->CreateContainer("cobjGFC", TList::Class(),AliAnalysisManager::kOutputContainer,outputGFC); 
607     mgr->ConnectInput(taskGFC,0,coutputFE); 
608     mgr->ConnectOutput(taskGFC,0,coutputGFC);
609     if (useWeights) {
610       mgr->ConnectInput(taskGFC,1,cinputWeights);
611       cinputWeights->SetData(weightsList);
612     } 
613   }
614   if(QC) {
615     TString outputQC = AliAnalysisManager::GetCommonFileName();
616     outputQC += ":outputQCanalysis";
617     outputQC+= type;
618
619     AliAnalysisDataContainer *coutputQC = mgr->CreateContainer("cobjQC", TList::Class(),AliAnalysisManager::kOutputContainer,outputQC); 
620     mgr->ConnectInput(taskQC,0,coutputFE); 
621     mgr->ConnectOutput(taskQC,0,coutputQC);
622     if (useWeights) {
623       mgr->ConnectInput(taskQC,1,cinputWeights);
624       cinputWeights->SetData(weightsList);
625     } 
626   }
627   if(FQD) {
628     TString outputFQD = AliAnalysisManager::GetCommonFileName();
629     outputFQD += ":outputFQDanalysis";
630     outputFQD+= type;
631     
632     AliAnalysisDataContainer *coutputFQD = mgr->CreateContainer("cobjFQD", TList::Class(),AliAnalysisManager::kOutputContainer,outputFQD); 
633     mgr->ConnectInput(taskFQD,0,coutputFE); 
634     mgr->ConnectOutput(taskFQD,0,coutputFQD);
635     if(useWeights) {
636       mgr->ConnectInput(taskFQD,1,cinputWeights);
637       cinputWeights->SetData(weightsList);
638     } 
639   }
640   if(MCEP) {
641     TString outputMCEP = AliAnalysisManager::GetCommonFileName();
642     outputMCEP += ":outputMCEPanalysis";
643     outputMCEP+= type;
644     
645     AliAnalysisDataContainer *coutputMCEP = mgr->CreateContainer("cobjMCEP", TList::Class(),AliAnalysisManager::kOutputContainer,outputMCEP); 
646     mgr->ConnectInput(taskMCEP,0,coutputFE); 
647     mgr->ConnectOutput(taskMCEP,0,coutputMCEP); 
648   }
649   
650
651   // Return analysis task
652   //===========================================================================
653   return taskFE;
654   
655
656
657 }
658
659
660
661
662