]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/macros/AODtrain.C
adding color reconnection option for Pythia6 generator (from Leticia Cunqueiro Mendez...
[u/mrichter/AliRoot.git] / ANALYSIS / macros / AODtrain.C
1 // ### Settings that make sense when using the Alien plugin
2 //==============================================================================
3 Int_t       runOnData          = 0;       // Set to 1 if processing real data
4 Int_t       iCollision         = 0;       // 0=pp, 1=Pb-Pb
5 Int_t       run_flag           = 1100;    // year (2011 = 1100)
6 //==============================================================================
7 Bool_t      doCDBconnect        =1;
8 Bool_t      usePhysicsSelection = kTRUE; // use physics selection
9 Bool_t      useTender           = kFALSE; // use tender wagon
10 Bool_t      useCentrality       = kTRUE; // centrality
11 Bool_t      useV0tender         = kFALSE;  // use V0 correction in tender
12 Bool_t      useDBG              = kTRUE;  // activate debugging
13 Bool_t      useMC               = kFALSE;  // use MC info
14 Bool_t      useKFILTER          = kFALSE;  // use Kinematics filter
15 Bool_t      useTR               = kFALSE;  // use track references
16 Bool_t      useCORRFW           = kFALSE; // do not change
17 Bool_t      useAODTAGS          = kFALSE; // use AOD tags
18 Bool_t      useSysInfo          = kFALSE; // use sys info
19
20 // ### Analysis modules to be included. Some may not be yet fully implemented.
21 //==============================================================================
22 Int_t       iAODhandler        = 1;      // Analysis produces an AOD or dAOD's
23 Int_t       iESDfilter         = 1;      // ESD to AOD filter (barrel + muon tracks)
24 Int_t       iMUONcopyAOD       = 1;      // Task that copies only muon events in a separate AOD (PWG3)
25 Int_t       iJETAN             = 0;      // Jet analysis (PWG4)
26 Int_t       iJETANdelta        = 0;      // Jet delta AODs
27 Int_t       iPWGHFvertexing     = 1;      // Vertexing HF task (PWG3)
28 Int_t       iPWGDQJPSIfilter    = 0;      // JPSI filtering (PWG3)
29 Int_t       iPWGHFd2h           = 1;      // D0->2 hadrons (PWG3)
30 Int_t       iPWGPP               =1;      // high pt filter task
31 Int_t       iPWGLFForward       = 1;      // Forward mult task (PWGLF)
32 Bool_t doPIDResponse  = 1;
33 Bool_t doPIDqa        = 1; //new
34
35 // ### Configuration macros used for each module
36 //==============================================================================
37  TString configPWGHFd2h = (iCollision==0)?"$ALICE_ROOT/PWGHF/vertexingHF/ConfigVertexingHF.C"
38                           :"$ALICE_ROOT/PWGHF/vertexingHF/ConfigVertexingHF_highmult.C";
39                                                  
40 // Temporaries.
41 void AODmerge();
42 void AddAnalysisTasks();
43 Bool_t LoadCommonLibraries();
44 Bool_t LoadAnalysisLibraries();
45 Bool_t LoadLibrary(const char *);
46 TChain *CreateChain();
47 const char *cdbPath = "raw://";
48 Int_t run_number = 0;
49
50 //______________________________________________________________________________
51 void AODtrain(Int_t merge=0)
52 {
53 // Main analysis train macro.
54
55   if (merge || doCDBconnect) {
56     TGrid::Connect("alien://");
57     if (!gGrid || !gGrid->IsConnected()) {
58       ::Error("QAtrain", "No grid connection");
59       return;
60     }
61   }
62   // Set temporary merging directory to current one
63   gSystem->Setenv("TMPDIR", gSystem->pwd());
64   // Set temporary compilation directory to current one
65   gSystem->SetBuildDir(gSystem->pwd(), kTRUE);
66    printf("==================================================================\n");
67    printf("===========    RUNNING FILTERING TRAIN   ==========\n");
68    printf("==================================================================\n");
69    printf("=  Configuring analysis train for:                               =\n");
70    if (usePhysicsSelection)   printf("=  Physics selection                                                =\n");
71    if (useTender)    printf("=  TENDER                                                        =\n");
72    if (iESDfilter)   printf("=  ESD filter                                                    =\n");
73    if (iMUONcopyAOD) printf("=  MUON copy AOD                                                 =\n");
74    if (iJETAN)       printf("=  Jet analysis                                                  =\n");
75    if (iJETANdelta)  printf("=     Jet delta AODs                                             =\n");
76    if (iPWGHFvertexing) printf("=  PWGHF vertexing                                                =\n");
77    if (iPWGDQJPSIfilter) printf("=  PWGDQ j/psi filter                                             =\n");
78    if (iPWGHFd2h) printf("=  PWGHF D0->2 hadrons QA                                     =\n");
79
80    // Load common libraries and set include path
81    if (!LoadCommonLibraries()) {
82       ::Error("AnalysisTrain", "Could not load common libraries");
83       return;
84    }
85     
86    // Make the analysis manager and connect event handlers
87    AliAnalysisManager *mgr  = new AliAnalysisManager("Analysis Train", "Production train");
88    if (useSysInfo) mgr->SetNSysInfo(100);
89    // Load analysis specific libraries
90    if (!LoadAnalysisLibraries()) {
91       ::Error("AnalysisTrain", "Could not load analysis libraries");
92       return;
93    }   
94
95    // Create input handler (input container created automatically)
96    // ESD input handler
97    AliESDInputHandler *esdHandler = new AliESDInputHandler();
98    mgr->SetInputEventHandler(esdHandler);       
99    // Monte Carlo handler
100    if (useMC) {
101       AliMCEventHandler* mcHandler = new AliMCEventHandler();
102       mgr->SetMCtruthEventHandler(mcHandler);
103       mcHandler->SetReadTR(useTR); 
104    }   
105    // AOD output container, created automatically when setting an AOD handler
106    if (iAODhandler) {
107       // AOD output handler
108       AliAODHandler* aodHandler   = new AliAODHandler();
109       aodHandler->SetOutputFileName("AliAOD.root");
110       mgr->SetOutputEventHandler(aodHandler);
111    }
112    // Debugging if needed
113    if (useDBG) mgr->SetDebugLevel(3);
114
115    AddAnalysisTasks(cdbPath);
116    if (merge) {
117       AODmerge();
118       mgr->InitAnalysis();
119       mgr->SetGridHandler(new AliAnalysisAlien);
120       mgr->StartAnalysis("gridterminate",0);
121       return;
122    }   
123    // Run the analysis                                                                                                                     
124    //
125    TChain *chain = CreateChain();
126    if (!chain) return;
127                                                                                                                                                    
128    TStopwatch timer;
129    timer.Start();
130    mgr->SetSkipTerminate(kTRUE);
131    if (mgr->InitAnalysis()) {
132       mgr->PrintStatus();
133       mgr->StartAnalysis("local", chain);
134    }
135    timer.Print();
136 }                                                                                                                                          
137                                                                                                                                             
138 //______________________________________________________________________________                                                           
139 void AddAnalysisTasks(const char *cdb_location){                                                                                                                                          
140   // Add all analysis task wagons to the train                                                                                               
141    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();                                                                     
142
143   //
144   // Tender and supplies. Needs to be called for every event.
145   //
146    AliAnalysisManager::SetCommonFileName("AODQA.root");
147    if (useTender) {
148       gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/TenderSupplies/AddTaskTender.C");
149       // IF V0 tender needed, put kTRUE below
150       AliAnalysisTaskSE *tender = AddTaskTender(useV0tender);
151 //      tender->SetDebugLevel(2);
152    }
153    //
154   // PIDResponse(JENS)
155   //
156   if (doPIDResponse) {
157     gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPIDResponse.C"); 
158     AliAnalysisTaskPIDResponse *PIDResponse = AddTaskPIDResponse();
159 //    PIDResponse->SelectCollisionCandidates(AliVEvent::kAny);
160   }  
161  
162   //
163   // PIDqa(JENS)
164   //
165   if (doPIDqa) {
166     gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPIDqa.C");
167     AliAnalysisTaskPIDqa *PIDQA = AddTaskPIDqa();
168     PIDQA->SelectCollisionCandidates(AliVEvent::kAny);
169   }  
170   // CDB connection
171   //
172   if (doCDBconnect && !useTender) {
173     gROOT->LoadMacro("$ALICE_ROOT/PWGPP/PilotTrain/AddTaskCDBconnect.C");
174     AliTaskCDBconnect *taskCDB = AddTaskCDBconnect(cdb_location, run_number);
175     if (!taskCDB) return;
176     AliCDBManager *cdb = AliCDBManager::Instance();
177     cdb->SetDefaultStorage(cdb_location);
178 //    taskCDB->SetRunNumber(run_number);
179   }    
180  
181    if (usePhysicsSelection) {
182    // Physics selection task
183       gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
184       mgr->RegisterExtraFile("event_stat.root");
185       AliPhysicsSelectionTask *physSelTask = AddTaskPhysicsSelection(useMC);
186 //      AliOADBPhysicsSelection * oadbDefaultPbPb = CreateOADBphysicsSelection();      
187 //      physSelTask->GetPhysicsSelection()->SetCustomOADBObjects(oadbDefaultPbPb,0,0);
188       mgr->AddStatisticsTask(AliVEvent::kAny);
189    }
190    
191
192 //Jacek
193    if (iPWGPP) {
194       gROOT->LoadMacro("$ALICE_ROOT/PWGPP/macros/AddTaskFilteredTree.C");
195       AddTaskFilteredTree("FilterEvents_Trees.root");
196    }   
197    
198    // Centrality (only Pb-Pb)
199    if (useCentrality) {
200       gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskCentrality.C");
201       AliCentralitySelectionTask *taskCentrality = AddTaskCentrality();
202       //taskCentrality->SelectCollisionCandidates(AliVEvent::kAny);
203    }
204    
205 // --- PWGLF - Forward (cholm@nbi.dk) -----------------------------
206    if (iPWGLFForward && usePhysicsSelection) { 
207         gROOT->LoadMacro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/AddTaskForwardMult.C");
208      UShort_t pwglfForwardSys = 0; // iCollision+1; // pp:1, PbPb:2, pPb:3
209      UShort_t pwglfSNN        = 0;            // GeV, 0==unknown
210      Short_t  pwglfField      = 0;
211      AddTaskForwardMult(useMC && useTR,        // Need track-refs 
212                         pwglfForwardSys,       // Collision system
213                         pwglfSNN, 
214                         pwglfField);
215         gROOT->LoadMacro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/AddTaskCentralMult.C");
216         AddTaskCentralMult(useMC, pwglfForwardSys, pwglfSNN, pwglfField);
217    }
218  
219     
220
221    if (iESDfilter) {
222       //  ESD filter task configuration.
223       gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/ESDfilter/macros/AddTaskESDFilter.C");
224       if (iMUONcopyAOD) {
225          printf("Registering delta AOD file\n");
226          mgr->RegisterExtraFile("AliAOD.Muons.root");
227           mgr->RegisterExtraFile("AliAOD.Dimuons.root");
228       }
229       AliAnalysisTaskESDfilter *taskesdfilter = 
230                  AddTaskESDFilter(useKFILTER, 
231                                   iMUONcopyAOD,         // write Muon AOD
232                                   kFALSE,               // write dimuon AOD 
233                                   kFALSE,               // usePhysicsSelection 
234                                   kFALSE,               // centrality OBSOLETE
235                                   kTRUE,                // enable TPS only tracks
236                                   kFALSE,               // disable cascades
237                                   kFALSE,               // disable kinks
238                                   run_flag);            // run flag (YY00)
239         mgr->RegisterExtraFile("AliAOD.Dimuons.root");
240    }   
241
242 // ********** PWG3 wagons ******************************************************           
243    // PWGHF vertexing
244    if (iPWGHFvertexing) {
245       gROOT->LoadMacro("$ALICE_ROOT/PWGHF/vertexingHF/macros/AddTaskVertexingHF.C");
246       if (!iPWGHFd2h) TFile::Cp(gSystem->ExpandPathName(configPWGHFd2h.Data()), "file:ConfigVertexingHF.C");
247       AliAnalysisTaskSEVertexingHF *taskvertexingHF = AddTaskVertexingHF();
248       if (!taskvertexingHF) ::Warning("AnalysisTrainNew", "AliAnalysisTaskSEVertexingHF cannot run for this train conditions - EXCLUDED");
249       else mgr->RegisterExtraFile("AliAOD.VertexingHF.root");
250       taskvertexingHF->SelectCollisionCandidates(0);
251    }   
252       
253    // PWGDQ JPSI filtering (only pp)
254    if (iPWGDQJPSIfilter && (iCollision==0)) {
255       gROOT->LoadMacro("$ALICE_ROOT/PWGDQ/dielectron/macros/AddTaskJPSIFilter.C");
256       AliAnalysisTaskSE *taskJPSIfilter = AddTaskJPSIFilter();
257       if (!taskJPSIfilter) ::Warning("AnalysisTrainNew", "AliAnalysisTaskDielectronFilter cannot run for this train conditions - EXCLUDED");
258       else mgr->RegisterExtraFile("AliAOD.Dielectron.root");
259       taskJPSIfilter->SelectCollisionCandidates(0);
260    }   
261
262    // PWGHF D2h
263    if (iPWGHFd2h) {   
264      gROOT->LoadMacro("$ALICE_ROOT/PWGHF/vertexingHF/AddD2HTrain.C");
265      TFile::Cp(gSystem->ExpandPathName(configPWGHFd2h.Data()), "file:ConfigVertexingHF.C");
266      AddD2HTrain(kFALSE, 1,0,0,0,0,0,0,0,0,0,0);                                 
267    }
268    
269    // ********** PWG4 wagons ******************************************************
270    // Jet analysis
271
272    // Configurations flags, move up?
273    TString kDeltaAODJetName = "AliAOD.Jets.root"; //
274    Bool_t  kIsPbPb = (iCollision==0)?false:true; // can be more intlligent checking the name of the data set
275    TString kDefaultJetBackgroundBranch = "";
276    TString kJetSubtractBranches = "";
277    UInt_t kHighPtFilterMask = 272;// from esd filter
278    UInt_t iPhysicsSelectionFlag = 0;
279    if (iJETAN) {
280      gROOT->LoadMacro("$ALICE_ROOT/PWGJE/macros/AddTaskJets.C");
281      // Default jet reconstructor running on ESD's
282      AliAnalysisTaskJets *taskjets = AddTaskJets("AOD","UA1",0.4,kHighPtFilterMask,1.,0); // no background subtraction     
283      if (!taskjets) ::Fatal("AnalysisTrainNew", "AliAnalysisTaskJets cannot run for this train conditions - EXCLUDED");
284      if(kDeltaAODJetName.Length()>0) taskjets->SetNonStdOutputFile(kDeltaAODJetName.Data());
285      if (iJETANdelta) {
286         //            AddTaskJetsDelta("AliAOD.Jets.root"); // need to modify this accordingly in the add task jets
287         mgr->RegisterExtraFile(kDeltaAODJetName.Data());
288         TString cTmp("");
289         if(kIsPbPb){
290           // UA1 intrinsic background subtraction
291           taskjets = AddTaskJets("AOD","UA1",0.4,kHighPtFilterMask,1.,2); // background subtraction
292           if(kDeltaAODJetName.Length()>0)taskjets->SetNonStdOutputFile(kDeltaAODJetName.Data());
293        }
294        // SICONE 
295        taskjets = AddTaskJets("AOD","SISCONE",0.4,kHighPtFilterMask,0.15,0); //no background subtration to be done later....                                                                                  
296        if(kDeltaAODJetName.Length()>0)taskjets->SetNonStdOutputFile(kDeltaAODJetName.Data());
297        cTmp = taskjets->GetNonStdBranch();
298        if(cTmp.Length()>0)kJetSubtractBranches += Form("%s ",cTmp.Data());
299          
300        // Add the clusters..
301        gROOT->LoadMacro("$ALICE_ROOT/PWGJE/macros/AddTaskJetCluster.C");
302        AliAnalysisTaskJetCluster *taskCl = 0;
303        Float_t fCenUp = 0;
304        Float_t fCenLo = 0;
305        Float_t fTrackEtaWindow = 0.9;
306        taskCl = AddTaskJetCluster("AOD","",kHighPtFilterMask,iPhysicsSelectionFlag,"KT",0.4,0,1, kDeltaAODJetName.Data(),0.15,fTrackEtaWindow,0); // this one is for the background and random jets, random cones with no skip                                                                                 
307        taskCl->SetBackgroundCalc(kTRUE);
308        taskCl->SetNRandomCones(10);
309        taskCl->SetCentralityCut(fCenLo,fCenUp);
310        taskCl->SetGhostEtamax(fTrackEtaWindow);
311        kDefaultJetBackgroundBranch = Form("%s_%s",AliAODJetEventBackground::StdBranchName(),taskCl->GetJetOutputBranch());
312
313        taskCl = AddTaskJetCluster("AOD","",kHighPtFilterMask,iPhysicsSelectionFlag,"ANTIKT",0.4,2,1,kDeltaAODJetName.Data(),0.15);
314        taskCl->SetCentralityCut(fCenLo,fCenUp);
315        if(kIsPbPb)taskCl->SetBackgroundBranch(kDefaultJetBackgroundBranch.Data());
316        taskCl->SetNRandomCones(10);
317        kJetSubtractBranches += Form("%s ",taskCl->GetJetOutputBranch());
318
319        taskCl = AddTaskJetCluster("AOD","",kHighPtFilterMask,iPhysicsSelectionFlag,"ANTIKT",0.2,0,1,kDeltaAODJetName.Data(),0.15);
320        taskCl->SetCentralityCut(fCenLo,fCenUp);
321        if(kIsPbPb)taskCl->SetBackgroundBranch(kDefaultJetBackgroundBranch.Data());
322        kJetSubtractBranches += Form("%s ",taskCl->GetJetOutputBranch());
323          
324        // DO THE BACKGROUND SUBTRACTION
325        if(kIsPbPb&&kJetSubtractBranches.Length()){
326          gROOT->LoadMacro("$ALICE_ROOT/PWGJE/macros/AddTaskJetBackgroundSubtract.C");
327          AliAnalysisTaskJetBackgroundSubtract *taskSubtract = 0;
328          taskSubtract = AddTaskJetBackgroundSubtract(kJetSubtractBranches,1,"B0","B%d");
329          taskSubtract->SetBackgroundBranch(kDefaultJetBackgroundBranch.Data());
330          if(kDeltaAODJetName.Length()>0)taskSubtract->SetNonStdOutputFile(kDeltaAODJetName.Data());
331        }
332      } 
333    }
334 }
335 //______________________________________________________________________________
336 Bool_t LoadCommonLibraries()
337 {
338 // Load common analysis libraries.
339    if (!gSystem->Getenv("ALICE_ROOT")) {
340       ::Error("AnalysisTrainNew.C::LoadCommonLibraries", "Analysis train requires that analysis libraries are compiled with a local AliRoot"); 
341       return kFALSE;
342    }   
343    Bool_t success = kTRUE;
344    // Load framework classes. Par option ignored here.
345    success &= LoadLibrary("libSTEERBase.so");
346    success &= LoadLibrary("libESD.so");
347    success &= LoadLibrary("libAOD.so");
348    success &= LoadLibrary("libANALYSIS.so");
349    success &= LoadLibrary("libOADB.so");
350    success &= LoadLibrary("libANALYSISalice.so");
351      success &= LoadLibrary("libESDfilter.so");
352    success &= LoadLibrary("libCORRFW.so");
353    gROOT->ProcessLine(".include $ALICE_ROOT/include");
354    if (success) {
355       ::Info("AnalysisTrainNew.C::LoadCommodLibraries", "Load common libraries:    SUCCESS");
356       ::Info("AnalysisTrainNew.C::LoadCommodLibraries", "Include path for Aclic compilation:\n%s",
357               gSystem->GetIncludePath());
358    } else {           
359       ::Info("AnalysisTrainNew.C::LoadCommodLibraries", "Load common libraries:    FAILED");
360    }   
361    return success;
362 }
363
364 //______________________________________________________________________________
365 Bool_t LoadAnalysisLibraries()
366 {
367 // Load common analysis libraries.
368    if (useTender || doCDBconnect) {
369       if (!LoadLibrary("TENDER") ||
370           !LoadLibrary("TENDERSupplies")) return kFALSE;
371    }       
372    // CDBconnect
373    if (doCDBconnect && !useTender) {
374       if (!LoadLibrary("PWGPP")) return kFALSE;
375    }
376           
377    if (iESDfilter || iPWGMuonTrain) {
378       if (!LoadLibrary("PWGmuon")) return kFALSE;
379    }   
380    // JETAN
381    if (iJETAN) {
382       if (!LoadLibrary("JETAN")) return kFALSE;
383    }
384    if (iJETANdelta) {
385       if (!LoadLibrary("JETAN") ||
386           !LoadLibrary("CGAL") ||
387           !LoadLibrary("fastjet") ||
388           !LoadLibrary("siscone") ||
389           !LoadLibrary("SISConePlugin") ||
390           !LoadLibrary("FASTJETAN")) return kFALSE;
391    }  
392    // PWG2 FORWARD
393    if (iPWGLFForward) {
394       //if (!LoadLibrary("PWGLFforward", mode, kTRUE)) return kFALSE;
395      if (!LoadLibrary("PWGLFforward2")) return kFALSE;
396   }   
397    
398    // PWG3 Vertexing HF
399    if (iPWGHFvertexing || iPWGHFd2h) {
400       if (!LoadLibrary("PWGflowBase") ||
401           !LoadLibrary("PWGflowTasks") ||
402           !LoadLibrary("PWGTRD") ||
403           !LoadLibrary("PWGHFvertexingHF")) return kFALSE;
404    }    
405  //    if (iPWGHFvertexing || iPWG3d2h) {
406  //     if (!LoadLibrary("PWG3base") ||
407  //         !LoadLibrary("PWGHFvertexingHF")) return kFALSE;
408  //  }   
409    // PWG3 dielectron
410    if (iPWGDQJPSIfilter) {
411       if (!LoadLibrary("PWGDQdielectron")) return kFALSE;
412    }   
413    
414    ::Info("AnalysisTrainNew.C::LoadAnalysisLibraries", "Load other libraries:   SUCCESS");
415    return kTRUE;
416 }
417
418 //______________________________________________________________________________
419 Bool_t LoadLibrary(const char *module)
420 {
421 // Load a module library in a given mode. Reports success.
422    Int_t result;
423    TString mod(module);
424    if (!mod.Length()) {
425       ::Error("AnalysisTrainNew.C::LoadLibrary", "Empty module name");
426       return kFALSE;
427    }   
428    // If a library is specified, just load it
429    if (mod.EndsWith(".so")) {
430       mod.Remove(mod.Index(".so"));
431       result = gSystem->Load(mod);
432       if (result < 0) {
433          ::Error("AnalysisTrainNew.C::LoadLibrary", "Could not load library %s", module);
434          return kFALSE;
435       }
436       return kTRUE;
437    } 
438    // Check if the library is already loaded
439    if (strlen(gSystem->GetLibraries(Form("%s.so", module), "", kFALSE)) > 0) return kTRUE;    
440    result = gSystem->Load(Form("lib%s.so", module));
441    if (result < 0) {
442       ::Error("AnalysisTrainNew.C::LoadLibrary", "Could not load module %s", module);
443       return kFALSE;
444    }
445    return kTRUE;
446 }           
447
448
449 //______________________________________________________________________________
450 TChain *CreateChain()
451 {
452 // Create the input chain
453    chain = new TChain("esdTree");
454    if (gSystem->AccessPathName("AliESDs.root")) 
455       ::Error("AnalysisTrainNew.C::CreateChain", "File: AliESDs.root not in ./data dir");
456    else 
457        chain->Add("AliESDs.root");
458    if (chain->GetNtrees()) return chain;
459    return NULL;
460 }   
461
462 //______________________________________________________________________________
463 void AODmerge()
464 {
465 // Merging method. No staging and no terminate phase.
466   TStopwatch timer;
467   timer.Start();
468   TString outputDir = "wn.xml";
469   TString outputFiles = "EventStat_temp.root,AODQA.root,AliAOD.root,AliAOD.VertexingHF.root,FilterEvents_Trees.root,AliAOD.Muons.root,AliAOD.Jets.root";
470   TString mergeExcludes = "";
471   TObjArray *list = outputFiles.Tokenize(",");
472   TIter *iter = new TIter(list);
473   TObjString *str;
474   TString outputFile;
475   Bool_t merged = kTRUE;
476   while((str=(TObjString*)iter->Next())) {
477     outputFile = str->GetString();
478     // Skip already merged outputs
479     if (!gSystem->AccessPathName(outputFile)) {
480        printf("Output file <%s> found. Not merging again.",outputFile.Data());
481        continue;
482     }
483     if (mergeExcludes.Contains(outputFile.Data())) continue;
484     merged = AliAnalysisAlien::MergeOutput(outputFile, outputDir, 10, 0);
485     if (!merged) {
486        printf("ERROR: Cannot merge %s\n", outputFile.Data());
487        continue;
488     }
489   }
490   // all outputs merged, validate
491   ofstream out;
492   out.open("outputs_valid", ios::out);
493   out.close();
494   timer.Print();
495 }