]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/macros/AddTaskFlowCentrality.C
ove all PID related things to AliFlowTrackCuts, make Bethe Bloch parameters user...
[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 // QA
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, kESD_TPConly, kESD_SPDtracklet
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 centralityName("");
80   centralityName+=Form("%.0f",centrMin);
81   centralityName+="-";
82   centralityName+=Form("%.0f",centrMax);
83
84   TString fileName(fileNameBase);
85   fileName.Append(".root");
86   //===========================================================================
87   printf("CREATE CUTS\n");
88   cout << "Used for RP: "<< rptypestr << endl;  
89   cout << "Used for POI: "<< poitypestr << endl;  
90   // EVENTS CUTS:
91   AliFlowEventCuts* cutsEvent = new AliFlowEventCuts("event cuts");
92   cutsEvent->SetCentralityPercentileRange(centrMin,centrMax);
93   cutsEvent->SetCentralityPercentileMethod(AliFlowEventCuts::kV0);
94   cutsEvent->SetRefMultMethod(AliFlowEventCuts::kV0);
95   //cutsEvent->SetCentralityPercentileMethod(AliFlowEventCuts::kSPD1tracklets);
96   cutsEvent->SetNContributorsRange(2);
97   cutsEvent->SetPrimaryVertexZrange(-10.,10.);
98   cutsEvent->SetCutSPDvertexerAnomaly(); //"Francesco's cut"
99   cutsEvent->SetCutZDCtiming();
100   
101   // RP TRACK CUTS:
102   AliFlowTrackCuts* cutsRP = new AliFlowTrackCuts("rp cuts");
103   cutsRP->SetParamType(rptype);
104   cutsRP->SetParamMix(rpmix);
105   cutsRP->SetPtRange(0.2,5.);
106   cutsRP->SetEtaRange(-0.8,0.8);
107   //cutsRP->SetRequireCharge(kTRUE);
108   //cutsRP->SetCharge(chargeRP);
109   //cutsRP->SetPID(PdgRP);
110   cutsRP->SetMinNClustersTPC(70);
111   cutsRP->SetMinChi2PerClusterTPC(0.1);
112   cutsRP->SetMaxChi2PerClusterTPC(4.0);
113   cutsRP->SetMinNClustersITS(2);
114   cutsRP->SetRequireITSRefit(kTRUE);
115   cutsRP->SetRequireTPCRefit(kTRUE);
116   //cutsRP->SetMaxChi2PerClusterITS(1.e+09);
117   cutsRP->SetMaxDCAToVertexXY(0.3);
118   cutsRP->SetMaxDCAToVertexZ(0.3);
119   //cutsRP->SetDCAToVertex2D(kTRUE);
120   //cutsRP->SetMaxNsigmaToVertex(1.e+10);
121   //cutsRP->SetRequireSigmaToVertex(kFALSE);
122   cutsRP->SetAcceptKinkDaughters(kFALSE);
123
124   // POI TRACK CUTS:
125   AliFlowTrackCuts* cutsPOI = new AliFlowTrackCuts("poi cuts");
126   cutsPOI->SetParamType(poitype);
127   cutsPOI->SetParamMix(poimix);
128   cutsPOI->SetPtRange(0.0,10.);
129   cutsPOI->SetEtaRange(-1.2,1.2);
130   //cutsPOI->SetRequireCharge(kTRUE);
131   //cutsPOI->SetCharge(chargeRP);
132   //cutsPOI->SetPID(PdgRP);
133   cutsPOI->SetMinNClustersTPC(80);
134   cutsPOI->SetMinChi2PerClusterTPC(0.1);
135   cutsPOI->SetMaxChi2PerClusterTPC(4.0);
136   cutsPOI->SetRequireITSRefit(kTRUE);
137   cutsPOI->SetRequireTPCRefit(kTRUE);
138   cutsPOI->SetMinNClustersITS(2);
139   //cutsPOI->SetMaxChi2PerClusterITS(1.e+09);
140   cutsPOI->SetMaxDCAToVertexXY(0.3);
141   cutsPOI->SetMaxDCAToVertexZ(0.3);
142   //cutsPOI->SetDCAToVertex2D(kTRUE);
143   //cutsPOI->SetMaxNsigmaToVertex(1.e+10);
144   //cutsPOI->SetRequireSigmaToVertex(kFALSE);
145   cutsPOI->SetAcceptKinkDaughters(kFALSE);
146   //cutsPOI->SetPID(AliPID::kProton, AliFlowTrackCuts::kTOFpid);
147   //cutsPOI->SetPID(AliPID::kPion, AliFlowTrackCuts::kTPCpid);
148   //cutsPOI->SetPID(AliPID::kProton, AliFlowTrackCuts::kTPCTOFpid);
149   //cutsPOI->SetTPCTOFpidCrossOverPt(0.4);
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   Bool_t useWeights  = WEIGHTS[0] || WEIGHTS[1] || WEIGHTS[2];
155   if (useWeights) cout<<"Weights are used"<<endl;
156   else cout<<"Weights are not used"<<endl;
157   
158   // Get the pointer to the existing analysis manager via the static access method.
159   //==============================================================================
160   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
161   if (!mgr) {
162     Error("AddTaskFlowEvent", "No analysis manager to connect to.");
163     return NULL;
164   }
165   
166   // Check the analysis type using the event handlers connected to the analysis
167   // manager. The availability of MC handler can also be checked here.
168   //==============================================================================
169   if (!mgr->GetInputEventHandler()) {
170     ::Error("AddTaskFlowEvent", "This task requires an input event handler");
171     return NULL;
172   }  
173
174   // Open external input files
175   //===========================================================================
176   //weights: 
177   TFile *weightsFile = NULL;
178   TList *weightsList = NULL;
179
180   if(useWeights) {
181     //open the file with the weights:
182     weightsFile = TFile::Open("weights.root","READ");
183     if(weightsFile) {
184       //access the list which holds the histos with weigths:
185       weightsList = (TList*)weightsFile->Get("weights");
186     }
187     else {
188       cout<<" WARNING: the file <weights.root> with weights from the previous run was not available."<<endl;
189       break;
190     } 
191   }
192   
193   //LYZ2
194   if (LYZ2SUM || LYZ2PROD) {
195     //read the outputfile of the first run
196     TString outputFileName = "AnalysisResults1.root";
197     TString pwd(gSystem->pwd());
198     pwd+="/";
199     pwd+=outputFileName.Data();
200     TFile *outputFile = NULL;
201     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {
202       cout<<"WARNING: You do not have an output file:"<<endl;
203       cout<<"         "<<pwd.Data()<<endl;
204       exit(0);
205     } else { outputFile = TFile::Open(pwd.Data(),"READ");}
206     
207     if (LYZ2SUM){  
208       // read the output directory from LYZ1SUM 
209       TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis" ;
210       inputFileNameLYZ2SUM += rptypestr;
211       cout<<"The input directory is "<<inputFileNameLYZ2SUM.Data()<<endl;
212       TFile* fInputFileLYZ2SUM = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2SUM.Data());
213       if(!fInputFileLYZ2SUM || fInputFileLYZ2SUM->IsZombie()) { 
214         cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ; 
215         break;
216       }
217       else {
218         TList* fInputListLYZ2SUM = (TList*)fInputFileLYZ2SUM->Get("cobjLYZ1SUM");
219         if (!fInputListLYZ2SUM) {cout<<"list is NULL pointer!"<<endl;}
220       }
221       cout<<"LYZ2SUM input file/list read..."<<endl;
222     }
223
224     if (LYZ2PROD){  
225       // read the output directory from LYZ1PROD 
226       TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis" ;
227       inputFileNameLYZ2PROD += rptypestr;
228       cout<<"The input directory is "<<inputFileNameLYZ2PROD.Data()<<endl;
229       TFile* fInputFileLYZ2PROD = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2PROD.Data());
230       if(!fInputFileLYZ2PROD || fInputFileLYZ2PROD->IsZombie()) { 
231         cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ; 
232         break;
233       }
234       else {
235         TList* fInputListLYZ2PROD = (TList*)fInputFileLYZ2PROD->Get("cobjLYZ1PROD");
236         if (!fInputListLYZ2PROD) {cout<<"list is NULL pointer!"<<endl;}
237       }
238       cout<<"LYZ2PROD input file/list read..."<<endl;
239     }
240   }
241
242   if (LYZEP) {
243     //read the outputfile of the second run
244     TString outputFileName = "AnalysisResults2.root";
245     TString pwd(gSystem->pwd());
246     pwd+="/";
247     pwd+=outputFileName.Data();
248     TFile *outputFile = NULL;
249     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {
250       cout<<"WARNING: You do not have an output file:"<<endl;
251       cout<<"         "<<pwd.Data()<<endl;
252       exit(0);
253     } else {
254       outputFile = TFile::Open(pwd.Data(),"READ");
255     }
256     
257     // read the output file from LYZ2SUM
258     TString inputFileNameLYZEP = "outputLYZ2SUManalysis" ;
259     inputFileNameLYZEP += rptypestr;
260     cout<<"The input file is "<<inputFileNameLYZEP.Data()<<endl;
261     TFile* fInputFileLYZEP = (TFile*)outputFile->FindObjectAny(inputFileNameLYZEP.Data());
262     if(!fInputFileLYZEP || fInputFileLYZEP->IsZombie()) { 
263       cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ; 
264       break;
265     }
266     else {
267       TList* fInputListLYZEP = (TList*)fInputFileLYZEP->Get("cobjLYZ2SUM");
268       if (!fInputListLYZEP) {cout<<"list is NULL pointer!"<<endl;}
269     }
270     cout<<"LYZEP input file/list read..."<<endl;
271   }
272   
273   
274   // Create the FMD task and add it to the manager
275   //===========================================================================
276   if (rptypestr == "FMD") {
277     AliFMDAnalysisTaskSE *taskfmd = NULL;
278     if (rptypestr == "FMD") {
279       taskfmd = new AliFMDAnalysisTaskSE("TaskFMD");
280       mgr->AddTask(taskfmd);
281       
282       AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
283       pars->Init();
284       pars->SetProcessPrimary(kTRUE); //for MC only
285       pars->SetProcessHits(kFALSE);
286       
287       //pars->SetRealData(kTRUE); //for real data
288       //pars->SetProcessPrimary(kFALSE); //for real data
289     }
290   }
291   
292   // Create the task, add it to the manager.
293   //===========================================================================
294   AliAnalysisTaskFlowEvent *taskFE = NULL;
295
296   if(useAfterBurner)
297     { 
298       taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent","",kFALSE,1);
299       taskFE->SetFlow(v1,v2,v3,v4); 
300       taskFE->SetNonFlowNumberOfTrackClones(numberOfTrackClones);
301       taskFE->SetAfterburnerOn();
302     }
303   else {taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent","",kFALSE); }
304   if (ExcludeRegion) {
305     taskFE->DefineDeadZone(excludeEtaMin, excludeEtaMax, excludePhiMin, excludePhiMax); 
306   }
307   taskFE->SetSubeventEtaRange(minA, maxA, minB, maxB);
308   if (UsePhysicsSelection) {
309     //taskFE->SelectCollisionCandidates(AliVEvent::kUserDefined);
310     taskFE->SelectCollisionCandidates(AliVEvent::kMB);
311     cout<<"Using Physics Selection"<<endl;
312   }
313   mgr->AddTask(taskFE);
314   
315   // Pass cuts for RPs and POIs to the task:
316   taskFE->SetCutsEvent(cutsEvent);
317
318   // Pass cuts for RPs and POIs to the task:
319   taskFE->SetCutsRP(cutsRP);
320   taskFE->SetCutsPOI(cutsPOI);
321
322   // Create the analysis tasks, add them to the manager.
323   //===========================================================================
324   if (SP){
325     AliAnalysisTaskScalarProduct *taskSP = new AliAnalysisTaskScalarProduct("TaskScalarProduct",WEIGHTS[0]);
326     taskSP->SetRelDiffMsub(1.0);
327     taskSP->SetApplyCorrectionForNUA(kTRUE);
328     mgr->AddTask(taskSP);
329   }
330   if (LYZ1SUM){
331     AliAnalysisTaskLeeYangZeros *taskLYZ1SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kTRUE);
332     taskLYZ1SUM->SetFirstRunLYZ(kTRUE);
333     taskLYZ1SUM->SetUseSumLYZ(kTRUE);
334     mgr->AddTask(taskLYZ1SUM);
335   }
336   if (LYZ1PROD){
337     AliAnalysisTaskLeeYangZeros *taskLYZ1PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kTRUE);
338     taskLYZ1PROD->SetFirstRunLYZ(kTRUE);
339     taskLYZ1PROD->SetUseSumLYZ(kFALSE);
340     mgr->AddTask(taskLYZ1PROD);
341   }
342   if (LYZ2SUM){
343     AliAnalysisTaskLeeYangZeros *taskLYZ2SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kFALSE);
344     taskLYZ2SUM->SetFirstRunLYZ(kFALSE);
345     taskLYZ2SUM->SetUseSumLYZ(kTRUE);
346     mgr->AddTask(taskLYZ2SUM);
347   }
348   if (LYZ2PROD){
349     AliAnalysisTaskLeeYangZeros *taskLYZ2PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kFALSE);
350     taskLYZ2PROD->SetFirstRunLYZ(kFALSE);
351     taskLYZ2PROD->SetUseSumLYZ(kFALSE);
352     mgr->AddTask(taskLYZ2PROD);
353   }
354   if (LYZEP){
355     AliAnalysisTaskLYZEventPlane *taskLYZEP = new AliAnalysisTaskLYZEventPlane("TaskLYZEventPlane");
356     mgr->AddTask(taskLYZEP);
357   }
358   if (GFC){
359     AliAnalysisTaskCumulants *taskGFC = new AliAnalysisTaskCumulants("TaskCumulants",useWeights);
360     taskGFC->SetUsePhiWeights(WEIGHTS[0]); 
361     taskGFC->SetUsePtWeights(WEIGHTS[1]);
362     taskGFC->SetUseEtaWeights(WEIGHTS[2]); 
363     mgr->AddTask(taskGFC);
364   }
365   if (QC){
366     AliAnalysisTaskQCumulants *taskQC = new AliAnalysisTaskQCumulants("TaskQCumulants",useWeights);
367     taskQC->SetUsePhiWeights(WEIGHTS[0]); 
368     taskQC->SetUsePtWeights(WEIGHTS[1]);
369     taskQC->SetUseEtaWeights(WEIGHTS[2]); 
370     taskQC->SetnBinsMult(10000);
371     taskQC->SetMinMult(0.);
372     taskQC->SetMaxMult(10000.);
373     taskQC->SetApplyCorrectionForNUA(kTRUE);
374     mgr->AddTask(taskQC);
375   }
376   if (FQD){
377     AliAnalysisTaskFittingQDistribution *taskFQD = new AliAnalysisTaskFittingQDistribution("TaskFittingQDistribution",kFALSE);
378     taskFQD->SetUsePhiWeights(WEIGHTS[0]); 
379     taskFQD->SetqMin(0.);
380     taskFQD->SetqMax(1000.);
381     taskFQD->SetqNbins(10000);
382     mgr->AddTask(taskFQD);
383   }
384   if (MCEP){
385     AliAnalysisTaskMCEventPlane *taskMCEP = new AliAnalysisTaskMCEventPlane("TaskMCEventPlane");
386     mgr->AddTask(taskMCEP);
387   }
388   if (MH){
389     AliAnalysisTaskMixedHarmonics *taskMH = new AliAnalysisTaskMixedHarmonics("TaskMixedHarmonics",useWeights);
390     taskMH->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]
391     taskMH->SetNoOfMultipicityBins(10000);
392     taskMH->SetMultipicityBinWidth(1.);
393     taskMH->SetMinMultiplicity(1.);
394     taskMH->SetCorrectForDetectorEffects(kTRUE);
395     taskMH->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>> (Remark: two nested loops)    
396     taskMH->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  
397     mgr->AddTask(taskMH);
398   }  
399   if (NL){
400     AliAnalysisTaskNestedLoops *taskNL = new AliAnalysisTaskNestedLoops("TaskNestedLoops",useWeights);
401     taskNL->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]
402     taskNL->SetEvaluateNestedLoopsForRAD(kTRUE); // RAD = Relative Angle Distribution
403     taskNL->SetEvaluateNestedLoopsForMH(kTRUE); // evalaute <<cos[n(phi1+phi2-2phi3)]>> (Remark: three nested loops)   
404     taskNL->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>>  (Remark: three nested loops)   
405     taskNL->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  
406     mgr->AddTask(taskNL);
407   }
408
409   // Create the output container for the data produced by the task
410   // Connect to the input and output containers
411   //===========================================================================
412   AliAnalysisDataContainer *cinput1 = mgr->GetCommonInputContainer();
413   
414   if (rptypestr == "FMD") {
415     AliAnalysisDataContainer *coutputFMD = 
416       mgr->CreateContainer(Form("BackgroundCorrected_%s",centralityName.Data()), TList::Class(), AliAnalysisManager::kExchangeContainer);
417     //input and output taskFMD     
418     mgr->ConnectInput(taskfmd, 0, cinput1);
419     mgr->ConnectOutput(taskfmd, 1, coutputFMD);
420     //input into taskFE
421     mgr->ConnectInput(taskFE,1,coutputFMD);
422   }
423   
424   AliAnalysisDataContainer *coutputFE = 
425   mgr->CreateContainer(Form("cobjFlowEventSimple_%s",centralityName.Data()),AliFlowEventSimple::Class(),AliAnalysisManager::kExchangeContainer);
426   mgr->ConnectInput(taskFE,0,cinput1); 
427   mgr->ConnectOutput(taskFE,1,coutputFE);
428   
429
430   // Create the output containers for the data produced by the analysis tasks
431   // Connect to the input and output containers
432   //===========================================================================
433   if (useWeights) {    
434     AliAnalysisDataContainer *cinputWeights = mgr->CreateContainer(Form("cobjWeights_%s",centralityName.Data()),
435                                                                    TList::Class(),AliAnalysisManager::kInputContainer); 
436   }
437
438   if(SP) {
439     TString outputSP = fileName;
440     outputSP += ":outputSPanalysis";
441     outputSP+= rptypestr;
442     AliAnalysisDataContainer *coutputSP = mgr->CreateContainer(Form("cobjSP_%s",centralityName.Data()), 
443                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputSP); 
444     mgr->ConnectInput(taskSP,0,coutputFE); 
445     mgr->ConnectOutput(taskSP,1,coutputSP); 
446     if (WEIGHTS[0]) {
447       mgr->ConnectInput(taskSP,1,cinputWeights);
448       cinputWeights->SetData(weightsList);
449     }
450   }
451   if(LYZ1SUM) {
452     TString outputLYZ1SUM = fileName;
453     outputLYZ1SUM += ":outputLYZ1SUManalysis";
454     outputLYZ1SUM+= rptypestr;
455     AliAnalysisDataContainer *coutputLYZ1SUM = mgr->CreateContainer(Form("cobjLYZ1SUM_%s",centralityName.Data()), 
456                                                                     TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1SUM); 
457     mgr->ConnectInput(taskLYZ1SUM,0,coutputFE);
458     mgr->ConnectOutput(taskLYZ1SUM,1,coutputLYZ1SUM);
459   }
460   if(LYZ1PROD) {
461     TString outputLYZ1PROD = fileName;
462     outputLYZ1PROD += ":outputLYZ1PRODanalysis";
463     outputLYZ1PROD+= rptypestr;
464     AliAnalysisDataContainer *coutputLYZ1PROD = mgr->CreateContainer(Form("cobjLYZ1PROD_%s",centralityName.Data()), 
465                                                                      TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1PROD); 
466     mgr->ConnectInput(taskLYZ1PROD,0,coutputFE); 
467     mgr->ConnectOutput(taskLYZ1PROD,1,coutputLYZ1PROD);
468   }
469   if(LYZ2SUM) {
470     AliAnalysisDataContainer *cinputLYZ2SUM = mgr->CreateContainer(Form("cobjLYZ2SUMin_%s",centralityName.Data()),
471                                                                    TList::Class(),AliAnalysisManager::kInputContainer);
472     TString outputLYZ2SUM = fileName;
473     outputLYZ2SUM += ":outputLYZ2SUManalysis";
474     outputLYZ2SUM+= rptypestr;
475     
476     AliAnalysisDataContainer *coutputLYZ2SUM = mgr->CreateContainer(Form("cobjLYZ2SUM_%s",centralityName.Data()), 
477                                                                     TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2SUM); 
478     mgr->ConnectInput(taskLYZ2SUM,0,coutputFE); 
479     mgr->ConnectInput(taskLYZ2SUM,1,cinputLYZ2SUM);
480     mgr->ConnectOutput(taskLYZ2SUM,1,coutputLYZ2SUM); 
481     cinputLYZ2SUM->SetData(fInputListLYZ2SUM);
482   }
483   if(LYZ2PROD) {
484     AliAnalysisDataContainer *cinputLYZ2PROD = mgr->CreateContainer(Form("cobjLYZ2PRODin_%s",centralityName.Data()),
485                                                                     TList::Class(),AliAnalysisManager::kInputContainer);
486     TString outputLYZ2PROD = fileName;
487     outputLYZ2PROD += ":outputLYZ2PRODanalysis";
488     outputLYZ2PROD+= rptypestr;
489     
490     AliAnalysisDataContainer *coutputLYZ2PROD = mgr->CreateContainer(Form("cobjLYZ2PROD_%s",centralityName.Data()), 
491                                                                      TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2PROD); 
492     mgr->ConnectInput(taskLYZ2PROD,0,coutputFE); 
493     mgr->ConnectInput(taskLYZ2PROD,1,cinputLYZ2PROD);
494     mgr->ConnectOutput(taskLYZ2PROD,1,coutputLYZ2PROD); 
495     cinputLYZ2PROD->SetData(fInputListLYZ2PROD);
496   }
497   if(LYZEP) {
498     AliAnalysisDataContainer *cinputLYZEP = mgr->CreateContainer(Form("cobjLYZEPin_%s",centralityName.Data()),
499                                                                  TList::Class(),AliAnalysisManager::kInputContainer);
500     TString outputLYZEP = fileName;
501     outputLYZEP += ":outputLYZEPanalysis";
502     outputLYZEP+= rptypestr;
503     
504     AliAnalysisDataContainer *coutputLYZEP = mgr->CreateContainer(Form("cobjLYZEP_%s",centralityName.Data()), 
505                                                                   TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZEP); 
506     mgr->ConnectInput(taskLYZEP,0,coutputFE); 
507     mgr->ConnectInput(taskLYZEP,1,cinputLYZEP);
508     mgr->ConnectOutput(taskLYZEP,1,coutputLYZEP); 
509     cinputLYZEP->SetData(fInputListLYZEP);
510   }
511   if(GFC) {
512     TString outputGFC = fileName;
513     outputGFC += ":outputGFCanalysis";
514     outputGFC+= rptypestr;
515     
516     AliAnalysisDataContainer *coutputGFC = mgr->CreateContainer(Form("cobjGFC_%s",centralityName.Data()), 
517                                                                 TList::Class(),AliAnalysisManager::kOutputContainer,outputGFC); 
518     mgr->ConnectInput(taskGFC,0,coutputFE); 
519     mgr->ConnectOutput(taskGFC,1,coutputGFC);
520     if (useWeights) {
521       mgr->ConnectInput(taskGFC,1,cinputWeights);
522       cinputWeights->SetData(weightsList);
523     } 
524   }
525   if(QC) {
526     TString outputQC = fileName;
527     outputQC += ":outputQCanalysis";
528     outputQC+= rptypestr;
529
530     AliAnalysisDataContainer *coutputQC = mgr->CreateContainer(Form("cobjQC_%s",centralityName.Data()), 
531                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputQC); 
532     mgr->ConnectInput(taskQC,0,coutputFE); 
533     mgr->ConnectOutput(taskQC,1,coutputQC);
534     if (useWeights) {
535       mgr->ConnectInput(taskQC,1,cinputWeights);
536       cinputWeights->SetData(weightsList);
537     }
538   }
539   if(FQD) {
540     TString outputFQD = fileName;
541     outputFQD += ":outputFQDanalysis";
542     outputFQD+= rptypestr;
543     
544     AliAnalysisDataContainer *coutputFQD = mgr->CreateContainer(Form("cobjFQD_%s",centralityName.Data()), 
545                                                                 TList::Class(),AliAnalysisManager::kOutputContainer,outputFQD); 
546     mgr->ConnectInput(taskFQD,0,coutputFE); 
547     mgr->ConnectOutput(taskFQD,1,coutputFQD);
548     if(useWeights) {
549       mgr->ConnectInput(taskFQD,1,cinputWeights);
550       cinputWeights->SetData(weightsList);
551     } 
552   }
553   if(MCEP) {
554     TString outputMCEP = fileName;
555     outputMCEP += ":outputMCEPanalysis";
556     outputMCEP+= rptypestr;
557     
558     AliAnalysisDataContainer *coutputMCEP = mgr->CreateContainer(Form("cobjMCEP_%s",centralityName.Data()), 
559                                                                  TList::Class(),AliAnalysisManager::kOutputContainer,outputMCEP); 
560     mgr->ConnectInput(taskMCEP,0,coutputFE);
561     mgr->ConnectOutput(taskMCEP,1,coutputMCEP); 
562   }
563   if(MH) {
564     TString outputMH = fileName;
565     outputMH += ":outputMHanalysis";
566     outputMH += rptypestr;
567         
568     AliAnalysisDataContainer *coutputMH = mgr->CreateContainer(Form("cobjMH_%s",centralityName.Data()), 
569                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputMH); 
570     mgr->ConnectInput(taskMH,0,coutputFE); 
571     mgr->ConnectOutput(taskMH,1,coutputMH); 
572     //if (useWeights) {
573     //  mgr->ConnectInput(taskMH,1,cinputWeights);
574     //  cinputWeights->SetData(weightsList);
575     //} 
576   }
577   if(NL) {
578     TString outputNL = fileName;
579     outputNL += ":outputNLanalysis";
580     outputNL += rptypestr;
581
582     AliAnalysisDataContainer *coutputNL = mgr->CreateContainer(Form("cobjNL_%s",centralityName.Data()), 
583                                                                TList::Class(),AliAnalysisManager::kOutputContainer,outputNL); 
584     mgr->ConnectInput(taskNL,0,coutputFE);
585     mgr->ConnectOutput(taskNL,1,coutputNL);
586     //if (useWeights) {
587     //  mgr->ConnectInput(taskNL,1,cinputWeights);
588     //  cinputWeights->SetData(weightsList);
589     //} 
590   }
591
592   ///////////////////////////////////////////////////////////////////////////////////////////
593   if (runQAtask)
594   {
595     AliAnalysisTaskQAflow* taskQAflow = new AliAnalysisTaskQAflow("TaskQAflow");
596     taskQAflow->SetEventCuts(cutsEvent);
597     taskQAflow->SetTrackCuts(cutsRP);
598     taskQAflow->SetFillNTuple(FillQAntuple);
599     taskQAflow->SetDoCorrelations(DoQAcorrelations);
600     mgr->AddTask(taskQAflow);
601     
602     Printf("centralityName %s",centralityName.Data());
603     TString taskQAoutputFileName(fileNameBase);
604     taskQAoutputFileName.Append("_QA.root");
605     AliAnalysisDataContainer* coutputQAtask = mgr->CreateContainer(Form("flowQA_%s",centralityName.Data()),
606                                               TObjArray::Class(),
607                                               AliAnalysisManager::kOutputContainer,
608                                               taskQAoutputFileName);
609     AliAnalysisDataContainer* coutputQAtaskTree = mgr->CreateContainer(Form("flowQAntuple_%s",centralityName.Data()),
610                                               TNtuple::Class(),
611                                               AliAnalysisManager::kOutputContainer,
612                                               taskQAoutputFileName);
613     mgr->ConnectInput(taskQAflow,0,mgr->GetCommonInputContainer());
614     mgr->ConnectInput(taskQAflow,1,coutputFE);
615     mgr->ConnectOutput(taskQAflow,1,coutputQAtask);
616     if (FillQAntuple) mgr->ConnectOutput(taskQAflow,2,coutputQAtaskTree);
617   }
618 }
619
620
621
622
623