3 * Example macro to run locally an analysis task for comparing the offline
4 * with the HLT esd tree.
6 * The output is a root file containing the histograms defined in the
7 * analysis task. There is one output file per detector.
9 * Run without arguments to get a few examples how to use the macro.
13 * aliroot -b -q -l compare_HLT_offline_local.C'("/home/blabla/AliESDs.root","global","./",kTRUE,10)' 2>&1 | tee task.log
14 * aliroot -b -q -l compare_HLT_offline_local.C'("/home/blabla/AliESDs.root","phos global pwg1",kTRUE,10)' 2>&1 | tee task.log
15 * aliroot -q compare-HLT-offline-local.C'("alien:///alice/data/2010/LHC10b/000115322/ESDs/pass1/10000115322040.20/AliESDs.root","global")' 2>&1 | tee log
18 * If alien:// is contained in the name of the file, then the macro connects to the grid to access the file.
20 * In case you want to run over many ESD files, then prepare a list of them in a .txt file and they will be chained for the analysis.
21 * The .txt file takes the place of the first argument in that case.
24 * @author Kalliopi.Kanaki@ift.uib.no, Hege.Erdal@student.uib.no
27 void compare_HLT_offline_local(TString file,
28 const char* detectorTask="global",
29 TString taskFolder="$ALICE_ROOT/HLT/QA/tasks/",
30 bool fUseHLTTrigger=kFALSE,
31 Long64_t nEvents=1234567890
38 gSystem->Load("libTree.so");
39 gSystem->Load("libGeom.so");
40 gSystem->Load("libVMC.so");
41 gSystem->Load("libPhysics.so");
43 //----------- Loading the required libraries ---------//
45 gSystem->Load("libSTEERBase.so");
46 gSystem->Load("libESD.so");
47 gSystem->Load("libAOD.so");
48 gSystem->Load("libANALYSIS.so");
49 gSystem->Load("libANALYSISalice.so");
50 gSystem->Load("libHLTbase.so");
52 gSystem->AddIncludePath("-I$ALICE_ROOT/PWG1/TPC -I.");
54 gSystem->Load("libTPCcalib.so");
55 gSystem->Load("libTRDbase.so");
56 gSystem->Load("libTRDrec.so");
57 gSystem->Load("libITSbase.so");
58 gSystem->Load("libITSrec.so");
59 gSystem->Load("libTENDER.so");
60 gSystem->Load("libPWG1.so");
62 gROOT->ProcessLine(".include $ALICE_ROOT/include");
63 //gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
65 Bool_t bPHOS=kFALSE, bGLOBAL=kFALSE, bEMCAL = kFALSE, bPWG1 = kFALSE, bD0=kFALSE;
67 TString allArgs = detectorTask;
70 TObjArray *pTokens = allArgs.Tokenize(" ");
72 for(int i=0; i<pTokens->GetEntries(); i++){
73 argument=((TObjString*)pTokens->At(i))->GetString();
74 if(argument.IsNull()) continue;
76 if(argument.CompareTo("phos", TString::kIgnoreCase)==0){
80 else if(argument.CompareTo("emcal", TString::kIgnoreCase)==0){
84 if(argument.CompareTo("global", TString::kIgnoreCase)==0){
88 if(argument.CompareTo("pwg1", TString::kIgnoreCase)==0){
92 if(argument.CompareTo("D0", TString::kIgnoreCase)==0){
96 if(argument.CompareTo("all",TString::kIgnoreCase)==0){
107 //-------------- Compile the analysis tasks ---------- //
110 gSystem->Load("libHLTbase");
111 gSystem->Load("libAliHLTUtil");
112 gSystem->Load("libAliHLTGlobal");
113 TString strTask1("AliAnalysisTaskCalo.cxx+");
114 TString strTask2("AliAnalysisTaskPHOS.cxx+");
115 gROOT->LoadMacro(taskFolder+strTask1);
116 gROOT->LoadMacro(taskFolder+strTask2);
117 cout << "\n========= You are loading the following tasks --> "<< taskFolder+strTask1 << " and " << taskFolder+strTask2 << endl;
121 gSystem->Load("libHLTbase");
122 gSystem->Load("libAliHLTUtil");
123 gSystem->Load("libAliHLTGlobal");
124 TString strTask1("AliAnalysisTaskCalo.cxx+");
125 TString strTask2("AliAnalysisTaskEMCAL.cxx+");
126 gROOT->LoadMacro(taskFolder+strTask1);
127 gROOT->LoadMacro(taskFolder+strTask2);
128 cout << "\n========= You are loading the following tasks --> "<< taskFolder+strTask1 << " and " << taskFolder+strTask2 << endl;
132 TString strTask("AliAnalysisTaskHLT.cxx+");
133 gROOT->LoadMacro(taskFolder+strTask);
134 cout << "\n========= You are loading the following task --> "<< taskFolder+strTask << endl;
137 TString strTask("AliAnalysisTaskD0Trigger.cxx+");
138 gROOT->LoadMacro(taskFolder+strTask);
139 cout << "\n========= You are loading the following task --> "<< taskFolder+strTask << endl;
142 if(bPWG1) gROOT->LoadMacro("$ALICE_ROOT/HLT/QA/tasks/macros/AddTaskPerformance.C");
144 if(file.Contains("alien")) TGrid::Connect("alien://");
146 if(file.Contains("AliESDs.root")){
147 TChain *chain = new TChain("esdTree");
151 //Constructs chain from filenames in *.txt
152 //on the form $DIR/AliESDs.root
153 else if(file.Contains(".txt")){
154 gROOT->LoadMacro("$ALICE_ROOT/PWG0/CreateESDChain.C");
155 chain=CreateESDChain(file.Data());
159 Printf("Chain is empty");
163 //To only select HLT triggered events
164 //Bool_t fUseHLTTrigger=kFALSE;
166 //-------- Make the analysis manager ---------------//
168 AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
169 AliESDInputHandler *esdH = new AliESDInputHandler;
171 //For the PWG1 task, setting HLT is handled inside AliPerformanceTask.C
172 if(!bPWG1) esdH->SetReadHLT();
173 esdH->SetReadFriends(kFALSE);
174 mgr->SetInputEventHandler(esdH);
175 mgr->SetNSysInfo(1000);
177 //AliPhysicsSelectionTask *physSelTask = AddTaskPhysicsSelection(kFALSE,kTRUE);
179 //-------------- define the tasks ------------//
182 AliAnalysisTaskHLTPHOS *taskPHOS = new AliAnalysisTaskHLTPHOS("offhlt_comparison_PHOS");
183 mgr->AddTask(taskPHOS);
184 AliAnalysisDataContainer *coutputPHOS = mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison.root");
185 mgr->ConnectInput(taskPHOS,0,mgr->GetCommonInputContainer());
186 mgr->ConnectOutput(taskPHOS,1,coutputPHOS);
190 AliAnalysisTaskHLTEMCAL *taskEMCAL = new AliAnalysisTaskHLTEMCAL("offhlt_comparison_EMCAL");
191 mgr->AddTask(taskEMCAL);
192 AliAnalysisDataContainer *coutputEMCAL = mgr->CreateContainer("emcal_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-EMCAL-comparison.root");
193 mgr->ConnectInput(taskEMCAL,0,mgr->GetCommonInputContainer());
194 mgr->ConnectOutput(taskEMCAL,1,coutputEMCAL);
198 AliAnalysisTaskHLT *taskGLOBAL = new AliAnalysisTaskHLT("offhlt_comparison_GLOBAL");
199 taskGLOBAL->SetUseHLTTriggerDecision(fUseHLTTrigger);
200 if(fUseHLTTrigger==kTRUE) printf("\n\nOnly HLT triggered events will be used to fill the distributions for task %s.\n\n", taskGLOBAL->GetName());
201 //taskGLOBAL->SelectCollisionCandidates();
202 mgr->AddTask(taskGLOBAL);
203 if(fUseHLTTrigger==kFALSE)
204 AliAnalysisDataContainer *coutputGLOBAL = mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison.root");
206 AliAnalysisDataContainer *coutputGLOBAL = mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer,"HLT-OFFLINE-GLOBAL-comparison_triggered.root");
207 mgr->ConnectInput(taskGLOBAL,0,mgr->GetCommonInputContainer());
208 mgr->ConnectOutput(taskGLOBAL,1,coutputGLOBAL);
213 // -- Add Task for HLT and Offline
214 AliPerformanceTask *HLTtpcQA = AddTaskPerformance(hasMC,kFALSE,kTRUE);
215 AliPerformanceTask *tpcQA = AddTaskPerformance(hasMC,kFALSE);
216 if(!HLTtpcQA || !tpcQA) {
217 Error("RunPerformanceTrain","AliPerformanceTask not created!");
222 float cuts[7]={0.5,0.04,0.7,0.8,0.05,-0.00025,0.7};
223 AliAnalysisTaskD0Trigger *taskD0 = new AliAnalysisTaskD0Trigger("offhlt_comparison_D0",cuts);
224 mgr->AddTask(taskD0);
225 AliAnalysisDataContainer *coutputD0 = mgr->CreateContainer("D0_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-D0-comparison.root");
226 mgr->ConnectInput(taskD0,0,mgr->GetCommonInputContainer());
227 mgr->ConnectOutput(taskD0,1,coutputD0);
230 if (!mgr->InitAnalysis()) return;
232 mgr->StartAnalysis("local",chain, nEvents);
238 void compare_HLT_offline_local(){
240 cout << " Usage examples:" << endl;
241 cout << " compare-HLT-offline-local.C'(file, taskOption, taskFolder, fUseHLTTrigger, nEvents)' 2>&1 | tee log" << endl;
242 cout << " compare-HLT-offline-local.C'(\"AliESDs.root\",\"global\")' 2>&1 | tee log" << endl;
243 cout << " compare-HLT-offline-local.C'(\"AliESDs.root\",\"global\",\"./\",kFALSE,nEvents)' 2>&1 | tee log" << endl;
244 cout << " compare-HLT-offline-local.C'(\"AliESDs.root\",\"global phos pwg1 D0\", \"./\", kTRUE, nEvents)' 2>&1 | tee log" << endl;
245 cout << " compare-HLT-offline-local.C'(\"files.txt\",\"all\")' 2>&1 | tee log" << endl;
246 cout << " compare-HLT-offline-local.C'(\"alien:///alice/data/2010/LHC10b/000115322/ESDs/pass1/10000115322040.20/AliESDs.root\",\"global\")' 2>&1 | tee log" << endl;