1 /////////////////////////////////////////////////////////////////////////////////////////////
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.
8 /////////////////////////////////////////////////////////////////////////////////////////////
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
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)
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;
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;
49 // For manipulating the event (for testing purposes)
50 const Bool_t AddToEvent = kFALSE;
51 Double_t ellipticflow = 0.05;
54 AliAnalysisTaskFlowEvent* AddTaskFlow(TString type, Bool_t* METHODS, Bool_t QA, Bool_t* WEIGHTS)
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];
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;
74 // Get the pointer to the existing analysis manager via the static access method.
75 //==============================================================================
76 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
78 Error("AddTaskFlowEvent", "No analysis manager to connect to.");
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");
90 // Open external input files
91 //===========================================================================
93 TFile *weightsFile = NULL;
94 TList *weightsList = NULL;
97 //open the file with the weights:
98 weightsFile = TFile::Open("weights.root","READ");
100 //access the list which holds the histos with weigths:
101 weightsList = (TList*)weightsFile->Get("weights");
104 cout<<" WARNING: the file <weights.root> with weights from the previous run was not available."<<endl;
109 if (LYZ2SUM || LYZ2PROD) {
110 //read the outputfile of the first run
111 TString outputFileName = "AnalysisResults1.root";
112 TString pwd(gSystem->pwd());
114 pwd+=outputFileName.Data();
115 TFile *outputFile = NULL;
116 if(gSystem->AccessPathName(pwd.Data(),kFileExists))
118 cout<<"WARNING: You do not have an output file:"<<endl;
119 cout<<" "<<pwd.Data()<<endl;
123 outputFile = TFile::Open(pwd.Data(),"READ");
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());
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 ;
139 TList* fInputListLYZ2SUM = (TList*)fInputFileLYZ2SUM->Get("cobjLYZ1SUM");
140 if (!fInputListLYZ2SUM) {cout<<"list is NULL pointer!"<<endl;}
142 cout<<"LYZ2SUM input file/list read..."<<endl;
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 ;
157 TList* fInputListLYZ2PROD = (TList*)fInputFileLYZ2PROD->Get("cobjLYZ1PROD");
158 if (!fInputListLYZ2PROD) {cout<<"list is NULL pointer!"<<endl;}
160 cout<<"LYZ2PROD input file/list read..."<<endl;
166 // read the output file from LYZ2SUM
167 TString outputFileName = "AnalysisResults2.root";
168 TString pwd(gSystem->pwd());
170 pwd+=outputFileName.Data();
171 TFile *outputFile = NULL;
172 if(gSystem->AccessPathName(pwd.Data(),kFileExists))
174 cout<<"WARNING: You do not have an output file:"<<endl;
175 cout<<" "<<pwd.Data()<<endl;
179 outputFile = TFile::Open(pwd.Data(),"READ");
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 ;
192 TList* fInputListLYZEP = (TList*)fInputFileLYZEP->Get("cobjLYZ2SUM");
193 if (!fInputListLYZEP) {cout<<"list is NULL pointer!"<<endl;}
195 cout<<"LYZEP input file/list read..."<<endl;
200 // Create the task, add it to the manager.
201 //===========================================================================
202 AliAnalysisTaskFlowEvent *taskFE = NULL;
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);
214 taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",kFALSE);
215 taskFE->SetAnalysisType(type);
216 taskFE->SetMinMult(multmin);
217 taskFE->SetMaxMult(multmax);
218 mgr->AddTask(taskFE);
221 // Create cuts using the correction framework (CORRFW)
222 //===========================================================================
224 //Set TList for the QA histograms
225 TList* qaRP = new TList();
226 TList* qaPOI = new TList();
229 //############# event cuts on multiplicity
230 AliCFEventGenCuts* mcEventCuts = new AliCFEventGenCuts("mcEventCuts","MC-level event cuts");
231 mcEventCuts->SetNTracksCut(multminESD,multmaxESD);
233 mcEventCuts->SetQAOn(qaRP);
235 AliCFEventRecCuts* recEventCuts = new AliCFEventRecCuts("recEventCuts","rec-level event cuts");
236 recEventCuts->SetNTracksCut(multminESD,multmaxESD);
238 recEventCuts->SetQAOn(qaRP);
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);
247 mcKineCuts1->SetQAOn(qaRP);
250 AliCFTrackKineCuts* mcKineCuts2 = new AliCFTrackKineCuts("mcKineCuts2","MC-level kinematic cuts");
251 mcKineCuts2->SetPtRange(ptmin2,ptmax2);
252 mcKineCuts2->SetRapidityRange(ymin2,ymax2);
253 //mcKineCuts2->SetChargeMC(charge2);
255 mcKineCuts2->SetQAOn(qaPOI);
258 AliCFParticleGenCuts* mcGenCuts1 = new AliCFParticleGenCuts("mcGenCuts1","MC particle generation cuts for RP");
259 mcGenCuts1->SetRequireIsPrimary();
260 if (UsePIDforRP) {mcGenCuts1->SetRequirePdgCode(PDG1);}
262 mcGenCuts1->SetQAOn(qaRP);
265 AliCFParticleGenCuts* mcGenCuts2 = new AliCFParticleGenCuts("mcGenCuts2","MC particle generation cuts for POI");
266 mcGenCuts2->SetRequireIsPrimary();
267 if (UsePIDforPOI) {mcGenCuts2->SetRequirePdgCode(PDG2);}
269 mcGenCuts2->SetQAOn(qaPOI);
272 //############# Acceptance Cuts
273 AliCFAcceptanceCuts *mcAccCuts1 = new AliCFAcceptanceCuts("mcAccCuts1","MC acceptance cuts");
274 mcAccCuts1->SetMinNHitITS(mintrackrefsITS1);
275 mcAccCuts1->SetMinNHitTPC(mintrackrefsTPC1);
277 mcAccCuts1->SetQAOn(qaRP);
280 AliCFAcceptanceCuts *mcAccCuts2 = new AliCFAcceptanceCuts("mcAccCuts2","MC acceptance cuts");
281 mcAccCuts2->SetMinNHitITS(mintrackrefsITS2);
282 mcAccCuts2->SetMinNHitTPC(mintrackrefsTPC2);
284 mcAccCuts2->SetQAOn(qaPOI);
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);
292 recKineCuts1->SetQAOn(qaRP);
295 AliCFTrackKineCuts *recKineCuts2 = new AliCFTrackKineCuts("recKineCuts2","rec-level kine cuts");
296 recKineCuts2->SetPtRange(ptmin2,ptmax2);
297 recKineCuts2->SetRapidityRange(ymin2,ymax2);
298 //recKineCuts2->SetChargeRec(charge2);
300 recKineCuts2->SetQAOn(qaPOI);
303 AliCFTrackQualityCuts *recQualityCuts1 = new AliCFTrackQualityCuts("recQualityCuts1","rec-level quality cuts");
304 recQualityCuts1->SetMinNClusterTPC(minclustersTPC1);
305 recQualityCuts1->SetStatus(AliESDtrack::kITSrefit);
307 recQualityCuts1->SetQAOn(qaRP);
309 AliCFTrackQualityCuts *recQualityCuts2 = new AliCFTrackQualityCuts("recQualityCuts2","rec-level quality cuts");
310 recQualityCuts2->SetMinNClusterTPC(minclustersTPC2);
311 recQualityCuts2->SetStatus(AliESDtrack::kITSrefit);
313 recQualityCuts2->SetQAOn(qaPOI);
316 AliCFTrackIsPrimaryCuts *recIsPrimaryCuts1 = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts1","rec-level isPrimary cuts");
317 recIsPrimaryCuts1->SetMaxNSigmaToVertex(maxnsigmatovertex1);
319 recIsPrimaryCuts1->SetQAOn(qaRP);
322 AliCFTrackIsPrimaryCuts *recIsPrimaryCuts2 = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts2","rec-level isPrimary cuts");
323 recIsPrimaryCuts2->SetMaxNSigmaToVertex(maxnsigmatovertex2);
325 recIsPrimaryCuts2->SetQAOn(qaPOI);
328 int n_species = AliPID::kSPECIES ;
329 Double_t* prior = new Double_t[n_species];
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 ;
337 AliCFTrackCutPid* cutPID1 = NULL;
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;
352 cutPID1->SetQAOn(qaRP);
356 AliCFTrackCutPid* cutPID2 = NULL;
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;
371 cutPID2->SetQAOn(qaPOI);
375 printf("CREATE EVENT CUTS\n");
376 TObjArray* mcEventList = new TObjArray(0);
377 mcEventList->AddLast(mcEventCuts);
379 TObjArray* recEventList = new TObjArray(0);
380 recEventList->AddLast(recEventCuts);
382 printf("CREATE MC KINE CUTS\n");
383 TObjArray* mcList1 = new TObjArray(0);
384 mcList1->AddLast(mcKineCuts1);
385 mcList1->AddLast(mcGenCuts1);
387 TObjArray* mcList2 = new TObjArray(0);
388 mcList2->AddLast(mcKineCuts2);
389 mcList2->AddLast(mcGenCuts2);
391 printf("CREATE ACCEPTANCE CUTS\n");
392 TObjArray* accList1 = new TObjArray(0) ;
393 accList1->AddLast(mcAccCuts1);
395 TObjArray* accList2 = new TObjArray(0) ;
396 accList2->AddLast(mcAccCuts2);
398 printf("CREATE RECONSTRUCTION CUTS\n");
399 TObjArray* recList1 = new TObjArray(0) ;
400 recList1->AddLast(recKineCuts1);
401 recList1->AddLast(recQualityCuts1);
402 recList1->AddLast(recIsPrimaryCuts1);
404 TObjArray* recList2 = new TObjArray(0) ;
405 recList2->AddLast(recKineCuts2);
406 recList2->AddLast(recQualityCuts2);
407 recList2->AddLast(recIsPrimaryCuts2);
409 printf("CREATE PID CUTS\n");
410 TObjArray* fPIDCutList1 = new TObjArray(0) ;
411 if(UsePIDforRP) {fPIDCutList1->AddLast(cutPID1);}
413 TObjArray* fPIDCutList2 = new TObjArray(0) ;
414 if (UsePIDforPOI) {fPIDCutList2->AddLast(cutPID2);}
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);
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);
438 taskFE->SetQAList1(qaRP);
439 taskFE->SetQAList2(qaPOI);
441 taskFE->SetCFManager1(cfmgr1);
442 taskFE->SetCFManager2(cfmgr2);
446 // Create the analysis tasks, add them to the manager.
447 //===========================================================================
449 AliAnalysisTaskScalarProduct *taskSP = new AliAnalysisTaskScalarProduct("TaskScalarProduct");
450 mgr->AddTask(taskSP);
453 AliAnalysisTaskLeeYangZeros *taskLYZ1SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kTRUE);
454 taskLYZ1SUM->SetFirstRunLYZ(kTRUE);
455 taskLYZ1SUM->SetUseSumLYZ(kTRUE);
456 mgr->AddTask(taskLYZ1SUM);
459 AliAnalysisTaskLeeYangZeros *taskLYZ1PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kTRUE);
460 taskLYZ1PROD->SetFirstRunLYZ(kTRUE);
461 taskLYZ1PROD->SetUseSumLYZ(kFALSE);
462 mgr->AddTask(taskLYZ1PROD);
465 AliAnalysisTaskLeeYangZeros *taskLYZ2SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kFALSE);
466 taskLYZ2SUM->SetFirstRunLYZ(kFALSE);
467 taskLYZ2SUM->SetUseSumLYZ(kTRUE);
468 mgr->AddTask(taskLYZ2SUM);
471 AliAnalysisTaskLeeYangZeros *taskLYZ2PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kFALSE);
472 taskLYZ2PROD->SetFirstRunLYZ(kFALSE);
473 taskLYZ2PROD->SetUseSumLYZ(kFALSE);
474 mgr->AddTask(taskLYZ2PROD);
477 AliAnalysisTaskLYZEventPlane *taskLYZEP = new AliAnalysisTaskLYZEventPlane("TaskLYZEventPlane");
478 mgr->AddTask(taskLYZEP);
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);
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);
495 AliAnalysisTaskFittingQDistribution *taskFQD = new AliAnalysisTaskFittingQDistribution("TaskFittingQDistribution",kFALSE);
496 taskFQD->SetUsePhiWeights(WEIGHTS[0]);
497 mgr->AddTask(taskFQD);
500 AliAnalysisTaskMCEventPlane *taskMCEP = new AliAnalysisTaskMCEventPlane("TaskMCEventPlane");
501 mgr->AddTask(taskMCEP);
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);
513 TString qaNameRPFE = AliAnalysisManager::GetCommonFileName();
514 qaNameRPFE += ":QAforRP_FE_";
517 AliAnalysisDataContainer *coutputQA1FE =
518 mgr->CreateContainer("QARPFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNameRPFE);
520 TString qaNamePOIFE = AliAnalysisManager::GetCommonFileName();
521 qaNamePOIFE += ":QAforPOI_FE_";
524 AliAnalysisDataContainer *coutputQA2FE =
525 mgr->CreateContainer("QAPOIFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNamePOIFE);
527 mgr->ConnectOutput(taskFE,1,coutputQA1FE);
528 mgr->ConnectOutput(taskFE,2,coutputQA2FE);
531 // Create the output containers for the data produced by the analysis tasks
532 // Connect to the input and output containers
533 //===========================================================================
535 AliAnalysisDataContainer *cinputWeights = mgr->CreateContainer("cobjWeights",TList::Class(),AliAnalysisManager::kInputContainer);
539 TString outputSP = AliAnalysisManager::GetCommonFileName();
540 outputSP += ":outputSPanalysis";
543 AliAnalysisDataContainer *coutputSP = mgr->CreateContainer("cobjSP", TList::Class(),AliAnalysisManager::kOutputContainer,outputSP);
544 mgr->ConnectInput(taskSP,0,coutputFE);
545 mgr->ConnectOutput(taskSP,0,coutputSP);
548 TString outputLYZ1SUM = AliAnalysisManager::GetCommonFileName();
549 outputLYZ1SUM += ":outputLYZ1SUManalysis";
550 outputLYZ1SUM+= type;
552 AliAnalysisDataContainer *coutputLYZ1SUM = mgr->CreateContainer("cobjLYZ1SUM", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1SUM);
553 mgr->ConnectInput(taskLYZ1SUM,0,coutputFE);
554 mgr->ConnectOutput(taskLYZ1SUM,0,coutputLYZ1SUM);
557 TString outputLYZ1PROD = AliAnalysisManager::GetCommonFileName();
558 outputLYZ1PROD += ":outputLYZ1PRODanalysis";
559 outputLYZ1PROD+= type;
561 AliAnalysisDataContainer *coutputLYZ1PROD = mgr->CreateContainer("cobjLYZ1PROD", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1PROD);
562 mgr->ConnectInput(taskLYZ1PROD,0,coutputFE);
563 mgr->ConnectOutput(taskLYZ1PROD,0,coutputLYZ1PROD);
566 AliAnalysisDataContainer *cinputLYZ2SUM = mgr->CreateContainer("cobjLYZ2SUMin",TList::Class(),AliAnalysisManager::kInputContainer);
567 TString outputLYZ2SUM = AliAnalysisManager::GetCommonFileName();
568 outputLYZ2SUM += ":outputLYZ2SUManalysis";
569 outputLYZ2SUM+= type;
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);
578 AliAnalysisDataContainer *cinputLYZ2PROD = mgr->CreateContainer("cobjLYZ2PRODin",TList::Class(),AliAnalysisManager::kInputContainer);
579 TString outputLYZ2PROD = AliAnalysisManager::GetCommonFileName();
580 outputLYZ2PROD += ":outputLYZ2PRODanalysis";
581 outputLYZ2PROD+= type;
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);
590 AliAnalysisDataContainer *cinputLYZEP = mgr->CreateContainer("cobjLYZEPin",TList::Class(),AliAnalysisManager::kInputContainer);
591 TString outputLYZEP = AliAnalysisManager::GetCommonFileName();
592 outputLYZEP += ":outputLYZEPanalysis";
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);
602 TString outputGFC = AliAnalysisManager::GetCommonFileName();
603 outputGFC += ":outputGFCanalysis";
606 AliAnalysisDataContainer *coutputGFC = mgr->CreateContainer("cobjGFC", TList::Class(),AliAnalysisManager::kOutputContainer,outputGFC);
607 mgr->ConnectInput(taskGFC,0,coutputFE);
608 mgr->ConnectOutput(taskGFC,0,coutputGFC);
610 mgr->ConnectInput(taskGFC,1,cinputWeights);
611 cinputWeights->SetData(weightsList);
615 TString outputQC = AliAnalysisManager::GetCommonFileName();
616 outputQC += ":outputQCanalysis";
619 AliAnalysisDataContainer *coutputQC = mgr->CreateContainer("cobjQC", TList::Class(),AliAnalysisManager::kOutputContainer,outputQC);
620 mgr->ConnectInput(taskQC,0,coutputFE);
621 mgr->ConnectOutput(taskQC,0,coutputQC);
623 mgr->ConnectInput(taskQC,1,cinputWeights);
624 cinputWeights->SetData(weightsList);
628 TString outputFQD = AliAnalysisManager::GetCommonFileName();
629 outputFQD += ":outputFQDanalysis";
632 AliAnalysisDataContainer *coutputFQD = mgr->CreateContainer("cobjFQD", TList::Class(),AliAnalysisManager::kOutputContainer,outputFQD);
633 mgr->ConnectInput(taskFQD,0,coutputFE);
634 mgr->ConnectOutput(taskFQD,0,coutputFQD);
636 mgr->ConnectInput(taskFQD,1,cinputWeights);
637 cinputWeights->SetData(weightsList);
641 TString outputMCEP = AliAnalysisManager::GetCommonFileName();
642 outputMCEP += ":outputMCEPanalysis";
645 AliAnalysisDataContainer *coutputMCEP = mgr->CreateContainer("cobjMCEP", TList::Class(),AliAnalysisManager::kOutputContainer,outputMCEP);
646 mgr->ConnectInput(taskMCEP,0,coutputFE);
647 mgr->ConnectOutput(taskMCEP,0,coutputMCEP);
651 // Return analysis task
652 //===========================================================================