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