- removed unnecessary variables and change the default detector task to 'global'.
[u/mrichter/AliRoot.git] / HLT / QA / tasks / macros / compare-HLT-offline-local.C
1 // $Id$
2 /*
3  * Example macro to run locally an analysis task for comparing the offline
4  * with the HLT esd tree.
5  *
6  * The output is a root file containing the histograms defined in the
7  * analysis task. There is one output file per detector.
8  *
9  * Usage:
10  * <pre>
11  *   aliroot -b -q -l compare_HLT_offline_local.C'("phos")' 2>&1 | tee task.log
12  *   aliroot -b -q -l compare_HLT_offline_local.C'("phos tpc global")' 2>&1 | tee task.log
13  *   aliroot -b -q -l compare_HLT_offline_local.C 2>&1 | tee task.log
14  * </pre>
15  *
16  * If no argument is specified, ALL detector tasks are run.
17  *
18  * @ingroup alihlt_tpc
19  * @author zbyin@mail.ccnu.edu.cn, Kalliopi.Kanaki@ift.uib.no
20  */
21
22 void compare_HLT_offline_local(const char* detectorTask="global"){
23  
24   TStopwatch timer;
25   timer.Start();
26
27   gSystem->Load("libTree.so");
28   gSystem->Load("libGeom.so");
29   gSystem->Load("libVMC.so");
30   gSystem->Load("libPhysics.so");
31  
32   //----------- Loading the required libraries ---------//
33
34   gSystem->Load("libSTEERBase.so");
35   gSystem->Load("libESD.so");
36   gSystem->Load("libAOD.so");
37   gSystem->Load("libANALYSIS.so");
38   gSystem->Load("libANALYSISalice.so");
39   gSystem->Load("libHLTbase.so");
40   gROOT->ProcessLine(".include $ALICE_ROOT/include");
41
42   
43   Bool_t bTPC=kFALSE, bPHOS=kFALSE, bITS=kFALSE, bGLOBAL=kFALSE;
44  
45   TString allArgs = detectorTask;
46   TString argument;
47  
48   TObjArray *pTokens = allArgs.Tokenize(" ");
49   if(pTokens){
50      for(int i=0; i<pTokens->GetEntries(); i++){
51          argument=((TObjString*)pTokens->At(i))->GetString();
52          if(argument.IsNull()) continue;
53
54          if(argument.CompareTo("tpc", TString::kIgnoreCase)==0){
55             bTPC = kTRUE;
56             continue;
57          }        
58          if(argument.CompareTo("phos", TString::kIgnoreCase)==0){
59             bPHOS = kTRUE;
60             continue;
61          }         
62          if(argument.CompareTo("its", TString::kIgnoreCase)==0){
63             bITS = kTRUE;
64             continue;
65          }      
66          if(argument.CompareTo("global", TString::kIgnoreCase)==0){
67             bGLOBAL = kTRUE;
68             continue;
69          }        
70          if(argument.CompareTo("all",TString::kIgnoreCase)==0){
71             bTPC    = kTRUE;
72             bPHOS   = kTRUE;
73             bITS    = kTRUE;
74             bGLOBAL = kTRUE;
75             continue;
76          }
77          else break;
78     }
79   }
80     
81   
82   //-------------- Compile the analysis tasks ---------- //
83   if(bTPC)    gROOT->LoadMacro("AliAnalysisTaskHLTTPC.cxx+"); 
84   if(bPHOS)   gROOT->LoadMacro("AliAnalysisTaskHLTPHOS.cxx+"); 
85   if(bITS)    gROOT->LoadMacro("AliAnalysisTaskHLTITS.cxx+");
86   if(bGLOBAL) gROOT->LoadMacro("AliAnalysisTaskHLT.cxx+");
87   
88   //TGrid::Connect("alien");
89
90   TChain *chain = new TChain("esdTree"); 
91   //chain->Add("alien:///alice/data/2010/LHC10b/000115322/10000115322040.110/AliESDs.root");
92   
93   chain->Add("~/7TeV/115322/10000115322040.110/AliESDs.root");
94   //chain->Add("...");
95    
96   //-------- Make the analysis manager ---------------//
97  
98   AliAnalysisManager *mgr  = new AliAnalysisManager("TestManager");
99   AliESDInputHandler *esdH = new AliESDInputHandler;
100   esdH->SetReadHLT();
101   mgr->SetInputEventHandler(esdH);  
102   mgr->SetNSysInfo(1000);
103  
104   //-------------- define the tasks ------------//
105   
106   if(bTPC){ 
107      AliAnalysisTaskHLTTPC *taskTPC = new AliAnalysisTaskHLTTPC("offhlt_comparison_TPC");
108      mgr->AddTask(taskTPC);
109      AliAnalysisDataContainer *coutput1 =  mgr->CreateContainer("tpc_histograms", TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-TPC-comparison.root");  
110      mgr->ConnectInput(taskTPC,0,mgr->GetCommonInputContainer());
111      mgr->ConnectOutput(taskTPC,1,coutput1);
112   }
113
114   if(bPHOS){
115      AliAnalysisTaskHLTPHOS *taskPHOS = new AliAnalysisTaskHLTPHOS("offhlt_comparison_PHOS");
116      mgr->AddTask(taskPHOS);
117      AliAnalysisDataContainer *coutput2 =  mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison.root");  
118      mgr->ConnectInput(taskPHOS,0,mgr->GetCommonInputContainer());
119      mgr->ConnectOutput(taskPHOS,1,coutput2);
120   }
121   
122   if(bITS){
123      AliAnalysisTaskHLTITS *taskITS = new AliAnalysisTaskHLTITS("offhlt_comparison_ITS");
124      mgr->AddTask(taskITS);
125      AliAnalysisDataContainer *coutput3 =  mgr->CreateContainer("its_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-ITS-comparison.root");  
126      mgr->ConnectInput(taskITS,0,mgr->GetCommonInputContainer());
127      mgr->ConnectOutput(taskITS,1,coutput3);
128   }
129  
130   if(bGLOBAL){
131      AliAnalysisTaskHLT *taskGLOBAL = new AliAnalysisTaskHLT("offhlt_comparison_GLOBAL");
132      mgr->AddTask(taskGLOBAL);
133      AliAnalysisDataContainer *coutput4 =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison.root");  
134      mgr->ConnectInput(taskGLOBAL,0,mgr->GetCommonInputContainer());
135      mgr->ConnectOutput(taskGLOBAL,1,coutput4);
136   }
137   
138   if (!mgr->InitAnalysis()) return;
139   mgr->PrintStatus();
140   mgr->StartAnalysis("local",chain);
141
142   timer.Stop();
143   timer.Print();
144 }