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