]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/macros/AddTaskFlow.C
afterburner
[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 // Define the range for eta subevents (for SP method)
11 //-----(FMD 1.7 - 5.0)-----
12 //Double_t minA = -5.0;
13 //Double_t maxA = -1.7;
14 //Double_t minB = 1.7;
15 //Double_t maxB = 5.0;
16 //-----(Tracklets 0.9 - 2.0)-----
17 //Double_t minA = -2.0;
18 //Double_t maxA = -0.9;
19 //Double_t minB = 0.9;
20 //Double_t maxB = 2.0;
21 //-----(Global 0.5 - 0.9)-----
22 Double_t minA = -0.9;
23 Double_t maxA = -0.5;
24 Double_t minB = 0.5;
25 Double_t maxB = 0.9;
26
27 // AFTERBURNER
28 Bool_t useAfterBurner=kFALSE;
29 Double_t v1=0.0;
30 Double_t v2=0.0;
31 Double_t v3=0.0;
32 Double_t v4=0.0;
33 Int_t numberOfTrackClones=0; //non-flow
34
35 // use physics selection class
36 Bool_t  UsePhysicsSelection = kTRUE;
37
38 // SETTING THE CUTS
39
40 //----------Event selection----------
41 Bool_t UseMultCutforESD = kTRUE;
42 //Bool_t UseMultCutforESD = kFALSE;
43 const Int_t multminESD = 1;  //used for CORRFW cuts 
44 const Int_t multmaxESD = 1000000; //used for CORRFW cuts 
45
46 Bool_t requireVtxCuts = kTRUE;
47 const Double_t vertexXmin = -1.e99; 
48 const Double_t vertexXmax = 1.e99;
49 const Double_t vertexYmin = -1.e99;
50 const Double_t vertexYmax = 1.e99;
51 const Double_t vertexZmin = -1.e99; 
52 const Double_t vertexZmax = 1.e99; 
53 const Int_t vertexNContributorsmin = 1;
54 const Int_t vertexNContributorsmax = 10000;
55
56 //Bool_t UseMultCut = kFALSE;
57 Bool_t UseMultCut = kTRUE;
58 const Int_t multmin = 1;     //used for AliFlowEventSimple (to set the centrality)
59 const Int_t multmax = 10000;     //used for AliFlowEventSimple (to set the centrality)
60 //const Int_t multmin = 10;     //used for AliFlowEventSimple (to set the centrality)
61 //const Int_t multmax = 1000000;     //used for AliFlowEventSimple (to set the centrality)
62
63
64 //----------For RP selection----------
65 // Use Global tracks ("Global"), or SPD tracklets ("Tracklet") 
66 // or FMD hits ("FMD") for the RP selection
67 const TString rptype = "Global";
68 //const TString rptype = "Tracklet";
69 //const TString rptype = "FMD";
70
71 //KINEMATICS (on generated and reconstructed tracks)
72 Bool_t UseKineforRP =  kTRUE;
73 const Double_t ptminRP = 0.0;
74 const Double_t ptmaxRP = 10.0;
75 const Double_t etaminRP  = -0.9;
76 const Double_t etamaxRP  = 0.9;
77 const Int_t    chargeRP = 1;  //not used
78 const Bool_t   isChargedRP = kTRUE;
79
80 //PID (on generated and reconstructed tracks)
81 Bool_t UsePIDforRP = kFALSE;
82 const Int_t PdgRP = 211;
83
84 //TRACK QUALITY (on reconstructed tracks only)
85 //see /CORRFW/AliCFTrackQualityCuts class
86 Bool_t UseTrackQualityforRP =  kTRUE;
87 const Int_t    minClustersTpcRP = 80;           //default = -1; 
88 const Double_t maxChi2PerClusterTpcRP = 4.0;    //default = 1.e+09;
89 const UShort_t minDedxClusterTpcRP = 0;
90 const Int_t    minClustersItsRP = 2;            //panos
91 const Double_t maxChi2PerClusterItsRP = 1.e+09; 
92 const Int_t    minClustersTrdRP = -1;
93 const Int_t    minTrackletTrdRP = -1;
94 const Int_t    minTrackletTrdPidRP = -1;
95 const Double_t maxChi2PerClusterTrdRP = 1.e+09;
96 const ULong_t  statusRP = AliESDtrack::kTPCrefit;   //AliESDtrack::kTPCrefit &  AliESDtrack::kITSrefit 
97
98 //PRIMARY (on reconstructed tracks only)
99 //see /CORRFW/AliCFTrackIsPrimaryCuts class
100 Bool_t UsePrimariesforRP = kTRUE;
101 const Bool_t   spdVertexRP = kFALSE;
102 const Bool_t   tpcVertexRP = kFALSE;
103 const Float_t  minDcaToVertexXyRP = 0.;
104 const Float_t  minDcaToVertexZRP = 0.;
105 const Float_t  maxDcaToVertexXyRP = 2.4;         //default = 1.e+10;  //2.4;
106 const Float_t  maxDcaToVertexZRP = 3.2;          //default = 1.e+10;  //3.2;
107 const Bool_t   dcaToVertex2dRP = kFALSE;         //default = kFALSE;
108 const Bool_t   absDcaToVertexRP = kTRUE;         //default = kTRUE;
109 const Double_t minNSigmaToVertexRP = 0.;
110 const Double_t maxNSigmaToVertexRP = 1.e+10; //3.; //1.e+10
111 const Double_t maxSigmaDcaXySP = 1.e+10;
112 const Double_t maxSigmaDcaZSP = 1.e+10;
113 const Bool_t   requireSigmaToVertexSP = kFALSE;
114 const Bool_t   acceptKinkDaughtersSP = kFALSE;  //default = kTRUE;
115
116 //ACCEPTANCE (on generated tracks only : AliMCParticle)
117 //see /CORRFW/AliCFAcceptanceCuts class
118 Bool_t UseAcceptanceforRP =  kFALSE; 
119 const Int_t  minTrackrefsItsRP = 0;//3;
120 const Int_t  minTrackrefsTpcRP = 0;//2;
121 const Int_t  minTrackrefsTrdRP = 0; 
122 const Int_t  minTrackrefsTofRP = 0; 
123 const Int_t  minTrackrefsMuonRP = 0; 
124 //default for all is 0
125
126 //----------For POI selection----------
127 //KINEMATICS (on generated and reconstructed tracks)
128 Bool_t UseKineforPOI = kTRUE;
129 const Double_t ptminPOI = 0.0;
130 const Double_t ptmaxPOI = 10.0;
131 const Double_t etaminPOI  = -0.5;
132 const Double_t etamaxPOI  = 0.5;
133 const Int_t    chargePOI = 1;  //not used
134 const Bool_t   isChargedPOI = kTRUE;
135
136 //PID (on generated and reconstructed tracks)
137 Bool_t UsePIDforPOI = kFALSE;
138 const Int_t PdgPOI = 321;
139
140 //TRACK QUALITY (on reconstructed tracks only)
141 //see /CORRFW/AliCFTrackQualityCuts class
142 Bool_t UseTrackQualityforPOI = kTRUE;
143 const Int_t    minClustersTpcPOI = 80;
144 const Double_t maxChi2PerClusterTpcPOI = 4.0;    
145 const UShort_t minDedxClusterTpcPOI = 0;
146 const Int_t    minClustersItsPOI = 2;
147 const Double_t maxChi2PerClusterItsPOI = 1.e+09;
148 const Int_t    minClustersTrdPOI = -1;
149 const Int_t    minTrackletTrdPOI = -1;
150 const Int_t    minTrackletTrdPidPOI = -1;
151 const Double_t maxChi2PerClusterTrdPOI = 1.e+09;
152 const ULong_t  statusPOI = AliESDtrack::kTPCrefit;   
153
154 //PRIMARY (on reconstructed tracks only)
155 //see /CORRFW/AliCFTrackIsPrimaryCuts class
156 Bool_t UsePrimariesforPOI = kTRUE;
157 const Bool_t   spdVertexPOI = kFALSE;
158 const Bool_t   tpcVertexPOI = kFALSE;
159 const Float_t  minDcaToVertexXyPOI = 0.;
160 const Float_t  minDcaToVertexZPOI = 0.;
161 const Float_t  maxDcaToVertexXyPOI = 2.4;
162 const Float_t  maxDcaToVertexZPOI = 3.2;
163 const Bool_t   dcaToVertex2dPOI =  kFALSE;
164 const Bool_t   absDcaToVertexPOI = kTRUE;
165 const Double_t minNSigmaToVertexPOI = 0.;
166 const Double_t maxNSigmaToVertexPOI = 1.e+10;  
167 const Double_t maxSigmaDcaXyPOI = 1.e+10;
168 const Double_t maxSigmaDcaZPOI = 1.e+10;
169 const Bool_t   requireSigmaToVertexPOI = kFALSE;
170 const Bool_t   acceptKinkDaughtersPOI = kFALSE;
171
172 //ACCEPTANCE (on generated tracks only : AliMCParticle)
173 //see /CORRFW/AliCFAcceptanceCuts class
174 Bool_t UseAcceptanceforPOI = kFALSE;
175 const Int_t minTrackrefsItsPOI = 3;
176 const Int_t minTrackrefsTpcPOI = 2;
177 const Int_t minTrackrefsTrdPOI = 0; 
178 const Int_t minTrackrefsTofPOI = 0; 
179 const Int_t minTrackrefsMuonPOI = 0; 
180
181 AliAnalysisTaskFlowEvent* AddTaskFlow(TString type, Bool_t* METHODS, Bool_t QA, Bool_t* WEIGHTS)
182 {
183   //boleans for the methods
184   Bool_t SP       = METHODS[0];
185   Bool_t LYZ1SUM  = METHODS[1];
186   Bool_t LYZ1PROD = METHODS[2];
187   Bool_t LYZ2SUM  = METHODS[3];
188   Bool_t LYZ2PROD = METHODS[4];
189   Bool_t LYZEP    = METHODS[5];
190   Bool_t GFC      = METHODS[6];
191   Bool_t QC       = METHODS[7];
192   Bool_t FQD      = METHODS[8];
193   Bool_t MCEP     = METHODS[9];      
194   Bool_t MH       = METHODS[10];
195   Bool_t NL       = METHODS[11];  
196   //for using weights
197   Bool_t useWeights  = WEIGHTS[0] || WEIGHTS[1] || WEIGHTS[2];
198   if (useWeights) cout<<"Weights are used"<<endl;
199   else cout<<"Weights are not used"<<endl;
200
201
202   // Get the pointer to the existing analysis manager via the static access method.
203   //==============================================================================
204   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
205   if (!mgr) {
206     Error("AddTaskFlowEvent", "No analysis manager to connect to.");
207     return NULL;
208   }   
209   
210   // Check the analysis type using the event handlers connected to the analysis
211   // manager. The availability of MC handler cann also be checked here.
212   //==============================================================================
213   if (!mgr->GetInputEventHandler()) {
214     ::Error("AddTaskFlowEvent", "This task requires an input event handler");
215     return NULL;
216   }  
217     
218   // Open external input files
219   //===========================================================================
220   //weights: 
221   TFile *weightsFile = NULL;
222   TList *weightsList = NULL;
223
224   if(useWeights) {
225     //open the file with the weights:
226     weightsFile = TFile::Open("weights.root","READ");
227     if(weightsFile) {
228       //access the list which holds the histos with weigths:
229       weightsList = (TList*)weightsFile->Get("weights");
230     }
231     else {
232       cout<<" WARNING: the file <weights.root> with weights from the previous run was not available."<<endl;
233       break;
234     } 
235   }
236     
237   //LYZ2
238   if (LYZ2SUM || LYZ2PROD) {
239     //read the outputfile of the first run
240     TString outputFileName = "AnalysisResults1.root";
241     TString pwd(gSystem->pwd());
242     pwd+="/";
243     pwd+=outputFileName.Data();
244     TFile *outputFile = NULL;
245     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {
246       cout<<"WARNING: You do not have an output file:"<<endl;
247       cout<<"         "<<pwd.Data()<<endl;
248       exit(0);
249     } else {
250       outputFile = TFile::Open(pwd.Data(),"READ");
251     }
252     
253     if (LYZ2SUM){  
254       // read the output directory from LYZ1SUM 
255       TString inputFileNameLYZ2SUM = "outputLYZ1SUManalysis" ;
256       inputFileNameLYZ2SUM += type;
257       cout<<"The input directory is "<<inputFileNameLYZ2SUM.Data()<<endl;
258       TFile* fInputFileLYZ2SUM = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2SUM.Data());
259       if(!fInputFileLYZ2SUM || fInputFileLYZ2SUM->IsZombie()) { 
260         cerr << " ERROR: To run LYZ2SUM you need the output file from LYZ1SUM. This file is not there! Please run LYZ1SUM first." << endl ; 
261         break;
262       }
263       else {
264         TList* fInputListLYZ2SUM = (TList*)fInputFileLYZ2SUM->Get("cobjLYZ1SUM");
265         if (!fInputListLYZ2SUM) {cout<<"list is NULL pointer!"<<endl;}
266       }
267       cout<<"LYZ2SUM input file/list read..."<<endl;
268     }
269
270     if (LYZ2PROD){  
271       // read the output directory from LYZ1PROD 
272       TString inputFileNameLYZ2PROD = "outputLYZ1PRODanalysis" ;
273       inputFileNameLYZ2PROD += type;
274       cout<<"The input directory is "<<inputFileNameLYZ2PROD.Data()<<endl;
275       TFile* fInputFileLYZ2PROD = (TFile*)outputFile->FindObjectAny(inputFileNameLYZ2PROD.Data());
276       if(!fInputFileLYZ2PROD || fInputFileLYZ2PROD->IsZombie()) { 
277         cerr << " ERROR: To run LYZ2PROD you need the output file from LYZ1PROD. This file is not there! Please run LYZ1PROD first." << endl ; 
278         break;
279       }
280       else {
281         TList* fInputListLYZ2PROD = (TList*)fInputFileLYZ2PROD->Get("cobjLYZ1PROD");
282         if (!fInputListLYZ2PROD) {cout<<"list is NULL pointer!"<<endl;}
283       }
284       cout<<"LYZ2PROD input file/list read..."<<endl;
285     }
286   }
287
288
289   if (LYZEP) {
290     //read the outputfile of the second run
291     TString outputFileName = "AnalysisResults2.root";
292     TString pwd(gSystem->pwd());
293     pwd+="/";
294     pwd+=outputFileName.Data();
295     TFile *outputFile = NULL;
296     if(gSystem->AccessPathName(pwd.Data(),kFileExists)) {
297       cout<<"WARNING: You do not have an output file:"<<endl;
298       cout<<"         "<<pwd.Data()<<endl;
299       exit(0);
300     } else {
301       outputFile = TFile::Open(pwd.Data(),"READ");
302     }
303
304     // read the output file from LYZ2SUM
305     TString inputFileNameLYZEP = "outputLYZ2SUManalysis" ;
306     inputFileNameLYZEP += type;
307     cout<<"The input file is "<<inputFileNameLYZEP.Data()<<endl;
308     TFile* fInputFileLYZEP = (TFile*)outputFile->FindObjectAny(inputFileNameLYZEP.Data());
309     if(!fInputFileLYZEP || fInputFileLYZEP->IsZombie()) { 
310       cerr << " ERROR: To run LYZEP you need the output file from LYZ2SUM. This file is not there! Please run LYZ2SUM first." << endl ; 
311       break;
312     }
313     else {
314       TList* fInputListLYZEP = (TList*)fInputFileLYZEP->Get("cobjLYZ2SUM");
315       if (!fInputListLYZEP) {cout<<"list is NULL pointer!"<<endl;}
316     }
317     cout<<"LYZEP input file/list read..."<<endl;
318   }
319   
320   
321   // Create the FMD task and add it to the manager
322   //===========================================================================
323   if (rptype == "FMD") {
324     AliFMDAnalysisTaskSE *taskfmd = NULL;
325     if (rptype == "FMD") {
326       taskfmd = new AliFMDAnalysisTaskSE("TaskFMD");
327       mgr->AddTask(taskfmd);
328   
329       AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
330       pars->Init();
331       pars->SetProcessPrimary(kTRUE);
332       pars->SetProcessHits(kFALSE);
333     }
334   }
335   
336
337   // Create the task, add it to the manager.
338   //===========================================================================
339   AliAnalysisTaskFlowEvent *taskFE = NULL;
340   if (QA) { 
341     if(useAfterBurner)
342     { 
343       taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",rptype,kTRUE,1);
344       taskFE->SetFlow(v1,v2,v3,v4); 
345       taskFE->SetNonFlowNumberOfTrackClones(numberOfTrackClones);
346     }
347     else {taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",rptype,kTRUE); }
348     taskFE->SetAnalysisType(type);
349     taskFE->SetRPType(rptype); //only for ESD
350     if (UseMultCut) {
351       taskFE->SetMinMult(multmin);
352       taskFE->SetMaxMult(multmax);
353     }
354     taskFE->SetSubeventEtaRange(minA, maxA, minB, maxB);
355     if (UsePhysicsSelection) {
356       taskFE->SelectCollisionCandidates();
357       cout<<"Using Physics Selection"<<endl;
358     }
359     mgr->AddTask(taskFE);
360   }
361   else { 
362     if(useAfterBurner)
363     { 
364       taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",rptype,kFALSE,1);
365       taskFE->SetFlow(v1,v2,v3,v4); 
366       taskFE->SetNonFlowNumberOfTrackClones(numberOfTrackClones);
367     }
368     else {taskFE = new AliAnalysisTaskFlowEvent("TaskFlowEvent",rptype,kFALSE); }
369     taskFE->SetAnalysisType(type);
370     if (UseMultCut) {
371       taskFE->SetMinMult(multmin);
372       taskFE->SetMaxMult(multmax);
373     }
374     taskFE->SetSubeventEtaRange(minA, maxA, minB, maxB);
375     if (UsePhysicsSelection) {
376       taskFE->SelectCollisionCandidates();
377       cout<<"Using Physics Selection"<<endl;
378     }
379     mgr->AddTask(taskFE);
380   }
381  
382   // Create cuts using the correction framework (CORRFW)
383   //===========================================================================
384   if (QA){
385     //Set TList for the QA histograms
386     TList* qaRP  = new TList(); 
387     TList* qaPOI = new TList();
388   }
389
390   //----------Event cuts----------
391   AliCFEventGenCuts* mcEventCuts = new AliCFEventGenCuts("mcEventCuts","MC-level event cuts");
392   mcEventCuts->SetNTracksCut(multminESD,multmaxESD); 
393   mcEventCuts->SetRequireVtxCuts(requireVtxCuts);
394   mcEventCuts->SetVertexXCut(vertexXmin, vertexXmax);
395   mcEventCuts->SetVertexYCut(vertexYmin, vertexYmax);
396   mcEventCuts->SetVertexZCut(vertexZmin, vertexZmax);
397   if (QA) { 
398     mcEventCuts->SetQAOn(qaRP);
399   }
400   AliCFEventRecCuts* recEventCuts = new AliCFEventRecCuts("recEventCuts","rec-level event cuts");
401   recEventCuts->SetNTracksCut(multminESD,multmaxESD); 
402   recEventCuts->SetRequireVtxCuts(requireVtxCuts);
403   recEventCuts->SetVertexXCut(vertexXmin, vertexXmax);
404   recEventCuts->SetVertexYCut(vertexYmin, vertexYmax);
405   recEventCuts->SetVertexZCut(vertexZmin, vertexZmax);
406   recEventCuts->SetVertexNContributors(vertexNContributorsmin,vertexNContributorsmax);
407   if (QA) { 
408     recEventCuts->SetQAOn(qaRP);
409   }
410   
411   //----------Cuts for RP----------
412   //KINEMATICS (MC and reconstructed)
413   AliCFTrackKineCuts* mcKineCutsRP = new AliCFTrackKineCuts("mcKineCutsRP","MC-level kinematic cuts");
414   mcKineCutsRP->SetPtRange(ptminRP,ptmaxRP);
415   mcKineCutsRP->SetEtaRange(etaminRP,etamaxRP);
416   //mcKineCutsRP->SetChargeMC(chargeRP);
417   mcKineCutsRP->SetRequireIsCharged(isChargedRP);
418   if (QA) { 
419     mcKineCutsRP->SetQAOn(qaRP);
420   }
421
422   AliCFTrackKineCuts *recKineCutsRP = new AliCFTrackKineCuts("recKineCutsRP","rec-level kine cuts");
423   recKineCutsRP->SetPtRange(ptminRP,ptmaxRP);
424   recKineCutsRP->SetEtaRange(etaminRP,etamaxRP);
425   //recKineCutsRP->SetChargeRec(chargeRP);
426   recKineCutsRP->SetRequireIsCharged(isChargedRP);
427   if (QA) { 
428     recKineCutsRP->SetQAOn(qaRP);
429   }
430
431   //PID (MC and reconstructed)
432   AliCFParticleGenCuts* mcGenCutsRP = new AliCFParticleGenCuts("mcGenCutsRP","MC particle generation cuts for RP");
433   mcGenCutsRP->SetRequireIsPrimary();
434   if (UsePIDforRP) {mcGenCutsRP->SetRequirePdgCode(PdgRP);}
435   if (QA) { 
436     mcGenCutsRP->SetQAOn(qaRP);
437   }
438
439   int n_species = AliPID::kSPECIES ;
440   Double_t* prior = new Double_t[n_species];
441   
442   prior[0] = 0.0244519 ;
443   prior[1] = 0.0143988 ;
444   prior[2] = 0.805747  ;
445   prior[3] = 0.0928785 ;
446   prior[4] = 0.0625243 ;
447   
448   AliCFTrackCutPid* cutPidRP = NULL;
449   if(UsePIDforRP) {
450     cutPidRP = new AliCFTrackCutPid("cutPidRP","ESD_PID for RP") ;
451     cutPidRP->SetPriors(prior);
452     cutPidRP->SetProbabilityCut(0.0);
453     cutPidRP->SetDetectors("TPC TOF");
454     switch(TMath::Abs(PDG1)) {
455     case 11   : cutPidRP->SetParticleType(AliPID::kElectron, kTRUE); break;
456     case 13   : cutPidRP->SetParticleType(AliPID::kMuon    , kTRUE); break;
457     case 211  : cutPidRP->SetParticleType(AliPID::kPion    , kTRUE); break;
458     case 321  : cutPidRP->SetParticleType(AliPID::kKaon    , kTRUE); break;
459     case 2212 : cutPidRP->SetParticleType(AliPID::kProton  , kTRUE); break;
460     default   : printf("UNDEFINED PID\n"); break;
461     }
462     if (QA) { 
463       cutPidRP->SetQAOn(qaRP); 
464     }
465   }
466   
467   //TRACK QUALITY
468   AliCFTrackQualityCuts *recQualityCutsRP = new AliCFTrackQualityCuts("recQualityCutsRP","rec-level quality cuts");
469   recQualityCutsRP->SetMinNClusterTPC(minClustersTpcRP);
470   //recQualityCutsRP->SetMinFoundClusterTPC(minFoundClustersTpcRP); //only for internal TPC QA
471   recQualityCutsRP->SetMaxChi2PerClusterTPC(maxChi2PerClusterTpcRP);
472   recQualityCutsRP->SetMinNdEdxClusterTPC(minDedxClusterTpcRP);     //to reject secondaries
473
474   recQualityCutsRP->SetMinNClusterITS(minClustersItsRP);
475   recQualityCutsRP->SetMaxChi2PerClusterITS(maxChi2PerClusterItsRP);
476
477   recQualityCutsRP->SetMinNClusterTRD(minClustersTrdRP);
478   recQualityCutsRP->SetMinNTrackletTRD(minTrackletTrdRP);
479   recQualityCutsRP->SetMinNTrackletTRDpid(minTrackletTrdPidRP);
480   recQualityCutsRP->SetMaxChi2PerTrackletTRD(maxChi2PerClusterTrdRP);
481   recQualityCutsRP->SetStatus(statusRP);  
482   if (QA) { 
483     recQualityCutsRP->SetQAOn(qaRP);
484   }
485
486   /* 
487   //How to set this?
488   void SetMaxCovDiagonalElements(Float_t c1=1.e+09, Float_t c2=1.e+09, Float_t c3=1.e+09, Float_t c4=1.e+09, Float_t c5=1.e+09)
489     {fCovariance11Max=c1;fCovariance22Max=c2;fCovariance33Max=c3;fCovariance44Max=c4;fCovariance55Max=c5;}
490   */
491
492   //PRIMARIES
493   AliCFTrackIsPrimaryCuts *recIsPrimaryCutsRP = new AliCFTrackIsPrimaryCuts("recIsPrimaryCutsRP","rec-level isPrimary cuts");
494   recIsPrimaryCutsRP->UseSPDvertex(spdVertexRP);
495   recIsPrimaryCutsRP->UseTPCvertex(tpcVertexRP);
496   recIsPrimaryCutsRP->SetMinDCAToVertexXY(minDcaToVertexXyRP); 
497   recIsPrimaryCutsRP->SetMinDCAToVertexZ(minDcaToVertexZRP);
498   recIsPrimaryCutsRP->SetMaxDCAToVertexXY(maxDcaToVertexXyRP);
499   recIsPrimaryCutsRP->SetMaxDCAToVertexZ(maxDcaToVertexZRP); 
500   recIsPrimaryCutsRP->SetDCAToVertex2D(dcaToVertex2dRP);
501   recIsPrimaryCutsRP->SetAbsDCAToVertex(absDcaToVertexRP);
502   recIsPrimaryCutsRP->SetMinNSigmaToVertex(minNSigmaToVertexRP); 
503   recIsPrimaryCutsRP->SetMaxNSigmaToVertex(maxNSigmaToVertexRP); 
504   recIsPrimaryCutsRP->SetMaxSigmaDCAxy(maxSigmaDcaXySP);
505   recIsPrimaryCutsRP->SetMaxSigmaDCAz(maxSigmaDcaZSP);
506   recIsPrimaryCutsRP->SetRequireSigmaToVertex(requireSigmaToVertexSP);
507   recIsPrimaryCutsRP->SetAcceptKinkDaughters(acceptKinkDaughtersSP);
508   if (QA) { 
509     recIsPrimaryCutsRP->SetQAOn(qaRP);
510   }
511   
512   //ACCEPTANCE
513   AliCFAcceptanceCuts *mcAccCutsRP = new AliCFAcceptanceCuts("mcAccCutsRP","MC acceptance cuts");
514   mcAccCutsRP->SetMinNHitITS(minTrackrefsItsRP);
515   mcAccCutsRP->SetMinNHitTPC(minTrackrefsTpcRP);
516   mcAccCutsRP->SetMinNHitTRD(minTrackrefsTrdRP); 
517   mcAccCutsRP->SetMinNHitTOF(minTrackrefsTofRP);
518   mcAccCutsRP->SetMinNHitMUON(minTrackrefsMuonRP);
519   if (QA) { 
520     mcAccCutsRP->SetQAOn(qaRP);
521   }
522
523   
524   //----------Cuts for POI----------
525   //KINEMATICS (MC and reconstructed)
526   AliCFTrackKineCuts* mcKineCutsPOI = new AliCFTrackKineCuts("mcKineCutsPOI","MC-level kinematic cuts");
527   mcKineCutsPOI->SetPtRange(ptminPOI,ptmaxPOI);
528   mcKineCutsPOI->SetEtaRange(etaminPOI,etamaxPOI);
529   //mcKineCutsPOI->SetChargeMC(chargePOI);
530   mcKineCutsPOI->SetRequireIsCharged(isChargedPOI);
531   if (QA) { 
532     mcKineCutsPOI->SetQAOn(qaPOI);
533   }
534   
535   AliCFTrackKineCuts *recKineCutsPOI = new AliCFTrackKineCuts("recKineCutsPOI","rec-level kine cuts");
536   recKineCutsPOI->SetPtRange(ptminPOI,ptmaxPOI);
537   recKineCutsPOI->SetEtaRange(etaminPOI,etamaxPOI);
538   //recKineCutsPOI->SetChargeRec(chargePOI);
539   recKineCutsPOI->SetRequireIsCharged(isChargedPOI);
540   if (QA) { 
541     recKineCutsPOI->SetQAOn(qaPOI);
542   }
543   
544   //PID (MC and reconstructed)
545   AliCFParticleGenCuts* mcGenCutsPOI = new AliCFParticleGenCuts("mcGenCutsPOI","MC particle generation cuts for POI");
546   mcGenCutsPOI->SetRequireIsPrimary();
547   if (UsePIDforPOI) {mcGenCutsPOI->SetRequirePdgCode(PdgPOI);}
548   if (QA) { 
549     mcGenCutsPOI->SetQAOn(qaPOI);
550   }
551
552   AliCFTrackCutPid* cutPidPOI = NULL;
553   if (UsePIDforPOI) {
554     cutPidPOI = new AliCFTrackCutPid("cutPidPOI","ESD_PID for POI") ;
555     cutPidPOI->SetPriors(prior);
556     cutPidPOI->SetProbabilityCut(0.0);
557     cutPidPOI->SetDetectors("TPC TOF");
558     switch(TMath::Abs(PDG2)) {
559     case 11   : cutPidPOI->SetParticleType(AliPID::kElectron, kTRUE); break;
560     case 13   : cutPidPOI->SetParticleType(AliPID::kMuon    , kTRUE); break;
561     case 211  : cutPidPOI->SetParticleType(AliPID::kPion    , kTRUE); break;
562     case 321  : cutPidPOI->SetParticleType(AliPID::kKaon    , kTRUE); break;
563     case 2212 : cutPidPOI->SetParticleType(AliPID::kProton  , kTRUE); break;
564     default   : printf("UNDEFINED PID\n"); break;
565     }
566     if (QA) { 
567       cutPidPOI->SetQAOn(qaPOI);
568     }
569   }
570
571   //TRACK QUALITY
572   AliCFTrackQualityCuts *recQualityCutsPOI = new AliCFTrackQualityCuts("recQualityCutsPOI","rec-level quality cuts");
573   recQualityCutsPOI->SetMinNClusterTPC(minClustersTpcPOI);
574   //recQualityCutsPOI->SetMinFoundClusterTPC(minFoundClustersTpcPOI); //only for internal TPC QA
575   recQualityCutsPOI->SetMaxChi2PerClusterTPC(maxChi2PerClusterTpcPOI);
576   recQualityCutsPOI->SetMinNdEdxClusterTPC(minDedxClusterTpcPOI);     //to reject secondaries
577
578   recQualityCutsPOI->SetMinNClusterITS(minClustersItsPOI);
579   recQualityCutsPOI->SetMaxChi2PerClusterITS(maxChi2PerClusterItsPOI);
580
581   recQualityCutsPOI->SetMinNClusterTRD(minClustersTrdPOI);
582   recQualityCutsPOI->SetMinNTrackletTRD(minTrackletTrdPOI);
583   recQualityCutsPOI->SetMinNTrackletTRDpid(minTrackletTrdPidPOI);
584   recQualityCutsPOI->SetMaxChi2PerTrackletTRD(maxChi2PerClusterTrdPOI);
585   recQualityCutsPOI->SetStatus(statusPOI); 
586   if (QA) { 
587     recQualityCutsPOI->SetQAOn(qaPOI);
588   }
589
590   //PRIMARIES
591   AliCFTrackIsPrimaryCuts *recIsPrimaryCutsPOI = new AliCFTrackIsPrimaryCuts("recIsPrimaryCutsPOI","rec-level isPrimary cuts");
592   recIsPrimaryCutsPOI->UseSPDvertex(spdVertexPOI);
593   recIsPrimaryCutsPOI->UseTPCvertex(tpcVertexPOI);
594   recIsPrimaryCutsPOI->SetMinDCAToVertexXY(minDcaToVertexXyPOI); 
595   recIsPrimaryCutsPOI->SetMinDCAToVertexZ(minDcaToVertexZPOI);
596   recIsPrimaryCutsPOI->SetMaxDCAToVertexXY(maxDcaToVertexXyPOI);
597   recIsPrimaryCutsPOI->SetMaxDCAToVertexZ(maxDcaToVertexZPOI); 
598   recIsPrimaryCutsPOI->SetDCAToVertex2D(dcaToVertex2dPOI);
599   recIsPrimaryCutsPOI->SetAbsDCAToVertex(absDcaToVertexPOI);
600   recIsPrimaryCutsPOI->SetMinNSigmaToVertex(minNSigmaToVertexPOI); 
601   recIsPrimaryCutsPOI->SetMaxNSigmaToVertex(maxNSigmaToVertexPOI); 
602   recIsPrimaryCutsPOI->SetMaxSigmaDCAxy(maxSigmaDcaXyPOI);
603   recIsPrimaryCutsPOI->SetMaxSigmaDCAz(maxSigmaDcaZPOI);
604   recIsPrimaryCutsPOI->SetRequireSigmaToVertex(requireSigmaToVertexPOI);
605   recIsPrimaryCutsPOI->SetAcceptKinkDaughters(acceptKinkDaughtersPOI);
606   if (QA) { 
607     recIsPrimaryCutsPOI->SetQAOn(qaPOI);
608   }
609
610   //ACCEPTANCE
611   AliCFAcceptanceCuts *mcAccCutsPOI = new AliCFAcceptanceCuts("mcAccCutsPOI","MC acceptance cuts");
612   mcAccCutsPOI->SetMinNHitITS(minTrackrefsItsPOI);
613   mcAccCutsPOI->SetMinNHitTPC(minTrackrefsTpcPOI);
614   mcAccCutsPOI->SetMinNHitTRD(minTrackrefsTrdPOI); 
615   mcAccCutsPOI->SetMinNHitTOF(minTrackrefsTofPOI);
616   mcAccCutsPOI->SetMinNHitMUON(minTrackrefsMuonPOI);
617   if (QA) { 
618     mcAccCutsPOI->SetQAOn(qaPOI);
619   }
620
621      
622   
623   //----------Create Cut Lists----------
624   printf("CREATE EVENT CUTS\n");
625   TObjArray* mcEventList = new TObjArray(0);  
626   if (UseMultCutforESD) mcEventList->AddLast(mcEventCuts);//cut on mult and vertex
627     
628   TObjArray* recEventList = new TObjArray(0);
629   if (UseMultCutforESD) recEventList->AddLast(recEventCuts);//cut on mult and vertex
630
631   printf("CREATE MC KINE CUTS\n");
632   TObjArray* mcListRP = new TObjArray(0);
633   if (UseKineforRP) mcListRP->AddLast(mcKineCutsRP); //cut on pt/eta/phi
634   mcListRP->AddLast(mcGenCutsRP); //cut on primary and if (UsePIDforRP) MC PID
635   
636   TObjArray* mcListPOI = new TObjArray(0);
637   if (UseKineforPOI) mcListPOI->AddLast(mcKineCutsPOI); //cut on pt/eta/phi
638   mcListPOI->AddLast(mcGenCutsPOI); //cut on primary and if (UsePIDforPOI) MC PID
639   
640   printf("CREATE MC ACCEPTANCE CUTS\n");
641   TObjArray* accListRP = new TObjArray(0) ;
642   if (UseAcceptanceforRP) accListRP->AddLast(mcAccCutsRP); //cut on number of track references
643   
644   TObjArray* accListPOI = new TObjArray(0) ;
645   if (UseAcceptanceforPOI) accListPOI->AddLast(mcAccCutsPOI); //cut on number of track references
646   
647   printf("CREATE ESD RECONSTRUCTION CUTS\n");
648   TObjArray* recListRP = new TObjArray(0) ;
649   if (UseTrackQualityforRP) recListRP->AddLast(recQualityCutsRP);   //track quality
650   if (UsePrimariesforRP)    recListRP->AddLast(recIsPrimaryCutsRP); //cut if it is a primary
651   if (UseKineforRP)         recListRP->AddLast(recKineCutsRP);      //cut on pt/eta/phi  
652
653   TObjArray* recListPOI = new TObjArray(0) ;
654   if (UseTrackQualityforPOI) recListPOI->AddLast(recQualityCutsPOI);   //track quality
655   if (UsePrimariesforPOI)    recListPOI->AddLast(recIsPrimaryCutsPOI); //cut if it is a primary
656   if (UseKineforPOI)         recListPOI->AddLast(recKineCutsPOI);      //cut on pt/eta/phi
657
658   printf("CREATE ESD PID CUTS\n");
659   TObjArray* fPIDCutListRP = new TObjArray(0) ;
660   if(UsePIDforRP) {fPIDCutListRP->AddLast(cutPidRP);} //cut on ESD PID
661   
662   TObjArray* fPIDCutListPOI = new TObjArray(0) ;
663   if (UsePIDforPOI)  {fPIDCutListPOI->AddLast(cutPidPOI);} //cut on ESD PID
664   
665
666   //----------Add Cut Lists to the CF Manager----------
667   printf("CREATE INTERFACE AND CUTS\n");
668   AliCFManager* cfmgrRP = new AliCFManager();
669   cfmgrRP->SetNStepEvent(3);
670   cfmgrRP->SetEventCutsList(AliCFManager::kEvtGenCuts,mcEventList); 
671   cfmgrRP->SetEventCutsList(AliCFManager::kEvtRecCuts,recEventList); 
672   cfmgrRP->SetNStepParticle(4); 
673   cfmgrRP->SetParticleCutsList(AliCFManager::kPartGenCuts,mcListRP);
674   cfmgrRP->SetParticleCutsList(AliCFManager::kPartAccCuts,accListRP);
675   cfmgrRP->SetParticleCutsList(AliCFManager::kPartRecCuts,recListRP);
676   cfmgrRP->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutListRP);
677   
678   AliCFManager* cfmgrPOI = new AliCFManager();
679   cfmgrPOI->SetNStepEvent(3);
680   cfmgrPOI->SetEventCutsList(AliCFManager::kEvtGenCuts,mcEventList); 
681   cfmgrPOI->SetEventCutsList(AliCFManager::kEvtRecCuts,recEventList); 
682   cfmgrPOI->SetNStepParticle(4); 
683   cfmgrPOI->SetParticleCutsList(AliCFManager::kPartGenCuts,mcListPOI);
684   cfmgrPOI->SetParticleCutsList(AliCFManager::kPartAccCuts,accListPOI);
685   cfmgrPOI->SetParticleCutsList(AliCFManager::kPartRecCuts,recListPOI);
686   cfmgrPOI->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutListPOI);
687   
688   if (QA) {
689     taskFE->SetQAList1(qaRP);
690     taskFE->SetQAList2(qaPOI);
691   }
692   taskFE->SetCFManager1(cfmgrRP);
693   taskFE->SetCFManager2(cfmgrPOI);
694
695
696
697   // Create the analysis tasks, add them to the manager.
698   //===========================================================================
699   if (SP){
700     AliAnalysisTaskScalarProduct *taskSP = new AliAnalysisTaskScalarProduct("TaskScalarProduct",WEIGHTS[0]);
701     taskSP->SetRelDiffMsub(1.0);
702     taskSP->SetApplyCorrectionForNUA(kFALSE);
703     mgr->AddTask(taskSP);
704   }
705   if (LYZ1SUM){
706     AliAnalysisTaskLeeYangZeros *taskLYZ1SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kTRUE);
707     taskLYZ1SUM->SetFirstRunLYZ(kTRUE);
708     taskLYZ1SUM->SetUseSumLYZ(kTRUE);
709     mgr->AddTask(taskLYZ1SUM);
710   }
711   if (LYZ1PROD){
712     AliAnalysisTaskLeeYangZeros *taskLYZ1PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kTRUE);
713     taskLYZ1PROD->SetFirstRunLYZ(kTRUE);
714     taskLYZ1PROD->SetUseSumLYZ(kFALSE);
715     mgr->AddTask(taskLYZ1PROD);
716   }
717   if (LYZ2SUM){
718     AliAnalysisTaskLeeYangZeros *taskLYZ2SUM = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosSUM",kFALSE);
719     taskLYZ2SUM->SetFirstRunLYZ(kFALSE);
720     taskLYZ2SUM->SetUseSumLYZ(kTRUE);
721     mgr->AddTask(taskLYZ2SUM);
722   }
723   if (LYZ2PROD){
724     AliAnalysisTaskLeeYangZeros *taskLYZ2PROD = new AliAnalysisTaskLeeYangZeros("TaskLeeYangZerosPROD",kFALSE);
725     taskLYZ2PROD->SetFirstRunLYZ(kFALSE);
726     taskLYZ2PROD->SetUseSumLYZ(kFALSE);
727     mgr->AddTask(taskLYZ2PROD);
728   }
729   if (LYZEP){
730     AliAnalysisTaskLYZEventPlane *taskLYZEP = new AliAnalysisTaskLYZEventPlane("TaskLYZEventPlane");
731     mgr->AddTask(taskLYZEP);
732   }
733   if (GFC){
734     AliAnalysisTaskCumulants *taskGFC = new AliAnalysisTaskCumulants("TaskCumulants",useWeights);
735     taskGFC->SetUsePhiWeights(WEIGHTS[0]); 
736     taskGFC->SetUsePtWeights(WEIGHTS[1]);
737     taskGFC->SetUseEtaWeights(WEIGHTS[2]); 
738     mgr->AddTask(taskGFC);
739   }
740   if (QC){
741     AliAnalysisTaskQCumulants *taskQC = new AliAnalysisTaskQCumulants("TaskQCumulants",useWeights);
742     taskQC->SetUsePhiWeights(WEIGHTS[0]); 
743     taskQC->SetUsePtWeights(WEIGHTS[1]);
744     taskQC->SetUseEtaWeights(WEIGHTS[2]); 
745     taskQC->SetnBinsMult(10000);
746     taskQC->SetMinMult(0.);
747     taskQC->SetMaxMult(10000.);
748     mgr->AddTask(taskQC);
749   }
750   if (FQD){
751     AliAnalysisTaskFittingQDistribution *taskFQD = new AliAnalysisTaskFittingQDistribution("TaskFittingQDistribution",kFALSE);
752     taskFQD->SetUsePhiWeights(WEIGHTS[0]); 
753     taskFQD->SetqMin(0.);
754     taskFQD->SetqMax(1000.);
755     taskFQD->SetqNbins(10000);
756     mgr->AddTask(taskFQD);
757   }
758   if (MCEP){
759     AliAnalysisTaskMCEventPlane *taskMCEP = new AliAnalysisTaskMCEventPlane("TaskMCEventPlane");
760     mgr->AddTask(taskMCEP);
761   }
762   if (MH){
763     AliAnalysisTaskMixedHarmonics *taskMH = new AliAnalysisTaskMixedHarmonics("TaskMixedHarmonics",useWeights);
764     taskMH->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]
765     taskMH->SetNoOfMultipicityBins(10);
766     taskMH->SetMultipicityBinWidth(2);
767     taskMH->SetMinMultiplicity(3);
768     taskMH->SetCorrectForDetectorEffects(kTRUE);
769     taskMH->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>> (Remark: two nested loops)    
770     taskMH->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  
771     mgr->AddTask(taskMH);
772   }  
773   if (NL){
774     AliAnalysisTaskNestedLoops *taskNL = new AliAnalysisTaskNestedLoops("TaskNestedLoops",useWeights);
775     taskNL->SetHarmonic(1); // n in cos[n(phi1+phi2-2phi3)] and cos[n(psi1+psi2-2phi3)]
776     taskNL->SetEvaluateNestedLoopsForRAD(kTRUE); // RAD = Relative Angle Distribution
777     taskNL->SetEvaluateNestedLoopsForMH(kTRUE); // evalaute <<cos[n(phi1+phi2-2phi3)]>> (Remark: three nested loops)   
778     taskNL->SetEvaluateDifferential3pCorrelator(kFALSE); // evaluate <<cos[n(psi1+psi2-2phi3)]>>  (Remark: three nested loops)   
779     taskNL->SetOppositeChargesPOI(kFALSE); // POIs psi1 and psi2 in cos[n(psi1+psi2-2phi3)] will have opposite charges  
780     mgr->AddTask(taskNL);
781   }
782   
783   // Create the output container for the data produced by the task
784   // Connect to the input and output containers
785   //===========================================================================
786   AliAnalysisDataContainer *cinput1 = mgr->GetCommonInputContainer();
787   
788   if (rptype == "FMD") {
789     AliAnalysisDataContainer *coutputFMD = 
790       mgr->CreateContainer("BackgroundCorrected", TList::Class(), AliAnalysisManager::kExchangeContainer);                        
791     //input and output taskFMD     
792     mgr->ConnectInput(taskfmd, 0, cinput1);
793     mgr->ConnectOutput(taskfmd, 1, coutputFMD);
794     //input into taskFE
795     mgr->ConnectInput(taskFE,1,coutputFMD);
796   }
797   
798   AliAnalysisDataContainer *coutputFE = 
799     mgr->CreateContainer("cobjFlowEventSimple",  AliFlowEventSimple::Class(),AliAnalysisManager::kExchangeContainer);
800   mgr->ConnectInput(taskFE,0,cinput1); 
801   mgr->ConnectOutput(taskFE,1,coutputFE);
802
803   if (QA) { 
804     TString qaNameRPFE = AliAnalysisManager::GetCommonFileName();
805     qaNameRPFE += ":QAforRP_FE_";
806     qaNameRPFE += type;
807
808     AliAnalysisDataContainer *coutputQA1FE =
809       mgr->CreateContainer("QARPFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNameRPFE); 
810     
811     TString qaNamePOIFE = AliAnalysisManager::GetCommonFileName();
812     qaNamePOIFE += ":QAforPOI_FE_";
813     qaNamePOIFE += type;
814         
815     AliAnalysisDataContainer *coutputQA2FE =
816       mgr->CreateContainer("QAPOIFE", TList::Class(),AliAnalysisManager::kOutputContainer,qaNamePOIFE); 
817
818     mgr->ConnectOutput(taskFE,2,coutputQA1FE); 
819     mgr->ConnectOutput(taskFE,3,coutputQA2FE); 
820   }
821
822   // Create the output containers for the data produced by the analysis tasks
823   // Connect to the input and output containers
824   //===========================================================================
825   if (useWeights) {    
826     AliAnalysisDataContainer *cinputWeights = mgr->CreateContainer("cobjWeights",TList::Class(),AliAnalysisManager::kInputContainer); 
827   }
828
829   if(SP) {
830     TString outputSP = AliAnalysisManager::GetCommonFileName();
831     outputSP += ":outputSPanalysis";
832     outputSP+= type;
833     
834     AliAnalysisDataContainer *coutputSP = mgr->CreateContainer("cobjSP", TList::Class(),AliAnalysisManager::kOutputContainer,outputSP); 
835     mgr->ConnectInput(taskSP,0,coutputFE); 
836     mgr->ConnectOutput(taskSP,1,coutputSP); 
837     if (WEIGHTS[0]) {
838       mgr->ConnectInput(taskSP,1,cinputWeights);
839       cinputWeights->SetData(weightsList);
840     } 
841   }
842   if(LYZ1SUM) {
843     TString outputLYZ1SUM = AliAnalysisManager::GetCommonFileName();
844     outputLYZ1SUM += ":outputLYZ1SUManalysis";
845     outputLYZ1SUM+= type;
846     
847     AliAnalysisDataContainer *coutputLYZ1SUM = mgr->CreateContainer("cobjLYZ1SUM", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1SUM); 
848     mgr->ConnectInput(taskLYZ1SUM,0,coutputFE); 
849     mgr->ConnectOutput(taskLYZ1SUM,1,coutputLYZ1SUM); 
850   }
851   if(LYZ1PROD) {
852     TString outputLYZ1PROD = AliAnalysisManager::GetCommonFileName();
853     outputLYZ1PROD += ":outputLYZ1PRODanalysis";
854     outputLYZ1PROD+= type;
855     
856     AliAnalysisDataContainer *coutputLYZ1PROD = mgr->CreateContainer("cobjLYZ1PROD", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ1PROD); 
857     mgr->ConnectInput(taskLYZ1PROD,0,coutputFE); 
858     mgr->ConnectOutput(taskLYZ1PROD,1,coutputLYZ1PROD);
859   }
860   if(LYZ2SUM) {
861     AliAnalysisDataContainer *cinputLYZ2SUM = mgr->CreateContainer("cobjLYZ2SUMin",TList::Class(),AliAnalysisManager::kInputContainer);
862     TString outputLYZ2SUM = AliAnalysisManager::GetCommonFileName();
863     outputLYZ2SUM += ":outputLYZ2SUManalysis";
864     outputLYZ2SUM+= type;
865     
866     AliAnalysisDataContainer *coutputLYZ2SUM = mgr->CreateContainer("cobjLYZ2SUM", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2SUM); 
867     mgr->ConnectInput(taskLYZ2SUM,0,coutputFE); 
868     mgr->ConnectInput(taskLYZ2SUM,1,cinputLYZ2SUM);
869     mgr->ConnectOutput(taskLYZ2SUM,1,coutputLYZ2SUM); 
870     cinputLYZ2SUM->SetData(fInputListLYZ2SUM);
871   }
872   if(LYZ2PROD) {
873     AliAnalysisDataContainer *cinputLYZ2PROD = mgr->CreateContainer("cobjLYZ2PRODin",TList::Class(),AliAnalysisManager::kInputContainer);
874     TString outputLYZ2PROD = AliAnalysisManager::GetCommonFileName();
875     outputLYZ2PROD += ":outputLYZ2PRODanalysis";
876     outputLYZ2PROD+= type;
877     
878     AliAnalysisDataContainer *coutputLYZ2PROD = mgr->CreateContainer("cobjLYZ2PROD", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZ2PROD); 
879     mgr->ConnectInput(taskLYZ2PROD,0,coutputFE); 
880     mgr->ConnectInput(taskLYZ2PROD,1,cinputLYZ2PROD);
881     mgr->ConnectOutput(taskLYZ2PROD,1,coutputLYZ2PROD); 
882     cinputLYZ2PROD->SetData(fInputListLYZ2PROD);
883   }
884   if(LYZEP) {
885     AliAnalysisDataContainer *cinputLYZEP = mgr->CreateContainer("cobjLYZEPin",TList::Class(),AliAnalysisManager::kInputContainer);
886     TString outputLYZEP = AliAnalysisManager::GetCommonFileName();
887     outputLYZEP += ":outputLYZEPanalysis";
888     outputLYZEP+= type;
889     
890     AliAnalysisDataContainer *coutputLYZEP = mgr->CreateContainer("cobjLYZEP", TList::Class(),AliAnalysisManager::kOutputContainer,outputLYZEP); 
891     mgr->ConnectInput(taskLYZEP,0,coutputFE); 
892     mgr->ConnectInput(taskLYZEP,1,cinputLYZEP);
893     mgr->ConnectOutput(taskLYZEP,1,coutputLYZEP); 
894     cinputLYZEP->SetData(fInputListLYZEP);
895   }
896   if(GFC) {
897     TString outputGFC = AliAnalysisManager::GetCommonFileName();
898     outputGFC += ":outputGFCanalysis";
899     outputGFC+= type;
900     
901     AliAnalysisDataContainer *coutputGFC = mgr->CreateContainer("cobjGFC", TList::Class(),AliAnalysisManager::kOutputContainer,outputGFC); 
902     mgr->ConnectInput(taskGFC,0,coutputFE); 
903     mgr->ConnectOutput(taskGFC,1,coutputGFC);
904     if (useWeights) {
905       mgr->ConnectInput(taskGFC,1,cinputWeights);
906       cinputWeights->SetData(weightsList);
907     } 
908   }
909   if(QC) {
910     TString outputQC = AliAnalysisManager::GetCommonFileName();
911     outputQC += ":outputQCanalysis";
912     outputQC+= type;
913
914     AliAnalysisDataContainer *coutputQC = mgr->CreateContainer("cobjQC", TList::Class(),AliAnalysisManager::kOutputContainer,outputQC); 
915     mgr->ConnectInput(taskQC,0,coutputFE); 
916     mgr->ConnectOutput(taskQC,1,coutputQC);
917     if (useWeights) {
918       mgr->ConnectInput(taskQC,1,cinputWeights);
919       cinputWeights->SetData(weightsList);
920     }    
921   }
922   if(FQD) {
923     TString outputFQD = AliAnalysisManager::GetCommonFileName();
924     outputFQD += ":outputFQDanalysis";
925     outputFQD+= type;
926     
927     AliAnalysisDataContainer *coutputFQD = mgr->CreateContainer("cobjFQD", TList::Class(),AliAnalysisManager::kOutputContainer,outputFQD); 
928     mgr->ConnectInput(taskFQD,0,coutputFE); 
929     mgr->ConnectOutput(taskFQD,1,coutputFQD);
930     if(useWeights) {
931       mgr->ConnectInput(taskFQD,1,cinputWeights);
932       cinputWeights->SetData(weightsList);
933     } 
934   }
935   if(MCEP) {
936     TString outputMCEP = AliAnalysisManager::GetCommonFileName();
937     outputMCEP += ":outputMCEPanalysis";
938     outputMCEP+= type;
939     
940     AliAnalysisDataContainer *coutputMCEP = mgr->CreateContainer("cobjMCEP", TList::Class(),AliAnalysisManager::kOutputContainer,outputMCEP); 
941     mgr->ConnectInput(taskMCEP,0,coutputFE); 
942     mgr->ConnectOutput(taskMCEP,1,coutputMCEP); 
943   }
944   if(MH) {
945     TString outputMH = AliAnalysisManager::GetCommonFileName();
946     outputMH += ":outputMHanalysis";
947     outputMH += type;
948         
949     AliAnalysisDataContainer *coutputMH = mgr->CreateContainer("cobjMH", TList::Class(),AliAnalysisManager::kOutputContainer,outputMH); 
950     mgr->ConnectInput(taskMH,0,coutputFE); 
951     mgr->ConnectOutput(taskMH,1,coutputMH); 
952     //if (useWeights) {
953     //  mgr->ConnectInput(taskMH,1,cinputWeights);
954     //  cinputWeights->SetData(weightsList);
955     //} 
956   }
957   if(NL) {
958     TString outputNL = AliAnalysisManager::GetCommonFileName();
959     outputNL += ":outputNLanalysis";
960     outputNL += type;
961         
962     AliAnalysisDataContainer *coutputNL = mgr->CreateContainer("cobjNL", TList::Class(),AliAnalysisManager::kOutputContainer,outputNL); 
963     mgr->ConnectInput(taskNL,0,coutputFE); 
964     mgr->ConnectOutput(taskNL,1,coutputNL); 
965     //if (useWeights) {
966     //  mgr->ConnectInput(taskNL,1,cinputWeights);
967     //  cinputWeights->SetData(weightsList);
968     //} 
969   }
970
971   // Return analysis task
972   //===========================================================================
973   return taskFE;
974   
975
976
977 }
978
979
980
981
982