]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/tasks/macros/compare-HLT-offline-local.C
- added the new function call that sets the pass for the centrality task
[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  * Run without arguments to get a few examples how to use the macro.
10  * 
11  * The arguments are:
12  * - the input file or txt file containing a list of ESDs to be processed (CreateESDChain takes 20 files as a default argument)
13  * - the task you want to use
14  * - the path of the task location
15  * - the beam type, "p-p" or "Pb-Pb", this is relevant ONLY for the central barrel task at the moment and is 
16  *   used to select proper binning and axes ranges for the THnSparse objects that it fills
17  * - options to make the central barrel task more flexible and lightweight; you can select if you want to 
18  *   fill the THnSparse object with only event or track properties or only HLT data or only OFF
19  *   possible options are: event-off event-hlt track-off track-hlt, all are turned on by default
20  * - boolean variable for selecting events which contain an HLT trigger
21  * - number of events to be analyzed
22  *
23  * If alien:// is placed before the input filename, then the macro connects to the grid to access the file.
24  * 
25  * 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.
26  * The .txt file takes the place of the first argument in that case.
27  *
28  * @ingroup alihlt_qa
29  * @author Kalliopi.Kanaki@ift.uib.no, Hege.Erdal@student.uib.no
30  */
31
32 void compare_HLT_offline_local( TString file
33                                ,const char* detectorTask="global"
34                                ,TString taskFolder="$ALICE_ROOT/HLT/QA/tasks/"
35                                ,TString beamType="p-p"
36                                ,TString options="event-off event-hlt track-off track-hlt"
37                                ,bool fUseHLTTrigger=kFALSE
38                                ,Long64_t nEvents=1234567890
39                               )
40 {
41
42   TStopwatch timer;
43   timer.Start();
44
45   gSystem->Load("libTree.so");
46   gSystem->Load("libGeom.so");
47   gSystem->Load("libVMC.so");
48   gSystem->Load("libPhysics.so");
49  
50   //----------- Loading the required libraries ---------//
51
52   gSystem->Load("libSTEERBase.so");
53   gSystem->Load("libESD.so");
54   gSystem->Load("libAOD.so");
55   gSystem->Load("libANALYSIS.so");
56   gSystem->Load("libANALYSISalice.so");
57   gSystem->Load("libHLTbase.so");
58  
59   gSystem->AddIncludePath("-I$ALICE_ROOT/HLT/BASE -I$ALICE_ROOT/PWG1/TPC -I. -I$ALICE_ROOT/STEER -I$ALICE_ROOT/ANALYSIS");
60   
61   gSystem->Load("libTPCcalib.so");
62   gSystem->Load("libTRDbase.so");
63   gSystem->Load("libTRDrec.so");
64   gSystem->Load("libITSbase.so");
65   gSystem->Load("libITSrec.so");
66   gSystem->Load("libTENDER.so");
67   gSystem->Load("libPWG1.so");
68  
69   gROOT->ProcessLine(".include $ALICE_ROOT/include");
70   //gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
71     
72   Bool_t bPHOS = kFALSE, bGLOBAL = kFALSE, bEMCAL = kFALSE, bPWG1 = kFALSE, bD0 = kFALSE, bCB = kFALSE;
73  
74   TString allArgs = detectorTask;
75   TString argument;
76  
77   TObjArray *pTokens = allArgs.Tokenize(" ");
78   if(pTokens){
79     for(int i=0; i<pTokens->GetEntries(); i++){
80       argument=((TObjString*)pTokens->At(i))->GetString();
81       if(argument.IsNull()) continue;
82
83       if(argument.CompareTo("phos", TString::kIgnoreCase)==0){
84         bPHOS = kTRUE;
85         continue;
86       }         
87       else if(argument.CompareTo("emcal", TString::kIgnoreCase)==0){
88         bEMCAL = kTRUE;
89         continue;
90       }         
91       if(argument.CompareTo("global", TString::kIgnoreCase)==0){
92         bGLOBAL = kTRUE;
93         continue;
94       }      
95       if(argument.CompareTo("pwg1", TString::kIgnoreCase)==0){
96         bPWG1 = kTRUE;
97         continue;
98       }
99       if(argument.CompareTo("D0", TString::kIgnoreCase)==0){
100         bD0 = kTRUE;
101         continue;
102       }
103       if(argument.CompareTo("cb", TString::kIgnoreCase)==0){
104         bCB = kTRUE;
105         continue;
106       } 
107       else break;
108     }
109   }
110       
111   //-------------- Compile the analysis tasks ---------- //
112   
113   if(bPHOS){
114     //gSystem->Load("libHLTbase");
115     gSystem->Load("libAliHLTUtil");
116     gSystem->Load("libAliHLTGlobal");
117     TString strTask1("AliAnalysisTaskHLTCalo.cxx+");
118     TString strTask2("AliAnalysisTaskHLTPHOS.cxx+");    
119     gROOT->LoadMacro(taskFolder+strTask1); 
120     gROOT->LoadMacro(taskFolder+strTask2); 
121     cout << "\n========= You are loading the following tasks --> "<< (taskFolder+strTask1).Chop()  << " and " <<  (taskFolder+strTask2).Chop() << endl;
122   }
123   
124   if(bEMCAL){
125     //gSystem->Load("libHLTbase");
126     gSystem->Load("libAliHLTUtil");
127     gSystem->Load("libAliHLTGlobal");
128     TString strTask1("AliAnalysisTaskHLTCalo.cxx+");
129     TString strTask2("AliAnalysisTaskHLTEMCAL.cxx+");
130     gROOT->LoadMacro(taskFolder+strTask1); 
131     gROOT->LoadMacro(taskFolder+strTask2); 
132     cout << "\n========= You are loading the following tasks --> "<< (taskFolder+strTask1).Chop()  << " and " <<  (taskFolder+strTask2).Chop() << endl;
133   }  
134   
135   if(bGLOBAL){
136      TString strTask("AliAnalysisTaskHLT.cxx+");
137      gROOT->LoadMacro(taskFolder+strTask);
138      cout << "\n========= You are loading the following task --> "<< (taskFolder+strTask).Chop()  << endl;
139   }
140   if(bD0){
141      TString strTask("AliAnalysisTaskD0Trigger.cxx+");
142      gROOT->LoadMacro(taskFolder+strTask); 
143      cout << "\n========= You are loading the following task --> "<< (taskFolder+strTask).Chop()  << endl;
144   }
145  
146   if(bCB){
147     TString strTask("AliAnalysisTaskHLTCentralBarrel.cxx+");
148     gROOT->LoadMacro(taskFolder+strTask);
149     cout << "\n========= You are loading the following task --> "<< (taskFolder+strTask).Chop()  << endl;
150   }
151   
152   if(bPWG1) gROOT->LoadMacro("$ALICE_ROOT/HLT/QA/tasks/macros/AddTaskPerformance.C");
153    
154   if(file.Contains("alien")) TGrid::Connect("alien://");
155     
156   if(file.Contains("AliESDs.root")){
157     TChain *chain = new TChain("esdTree"); 
158     chain->Add(file);
159   }
160   
161   //Constructs chain from filenames in *.txt
162   //on the form $DIR/AliESDs.root
163   else if(file.Contains(".txt")){
164     gROOT->LoadMacro("$ALICE_ROOT/PWG0/CreateESDChain.C");
165     chain=CreateESDChain(file.Data(),200);
166   }
167   
168   else if(!file){
169     printf("File %s does not exist or is corrupted.\n",file.Data());
170     return;  
171   }
172
173   if(!chain){
174     Printf("Chain is empty.\n");
175     return;
176   }
177    
178   //-------- Make the analysis manager ---------------//
179  
180   AliAnalysisManager *mgr  = new AliAnalysisManager("TestManager");
181   AliESDInputHandler *esdH = new AliESDInputHandler;
182
183   //For the PWG1 task, setting HLT is handled inside AliPerformanceTask.C
184   if(!bPWG1)  esdH->SetReadHLT();
185   esdH->SetReadFriends(kFALSE);
186   mgr->SetInputEventHandler(esdH);  
187   mgr->SetNSysInfo(1000);
188
189   //AliPhysicsSelectionTask *physSelTask = AddTaskPhysicsSelection(kFALSE,kTRUE);
190  
191   //-------------- define the tasks ------------//
192   
193   if(bPHOS){
194     AliAnalysisTaskHLTPHOS *taskPHOS = new AliAnalysisTaskHLTPHOS("offhlt_comparison_PHOS");
195     taskPHOS->SetUseHLTTriggerDecision(fUseHLTTrigger);
196     if(fUseHLTTrigger==kTRUE) printf("\n\nOnly HLT triggered events will be used to fill the distributions for task %s.\n\n", taskPHOS->GetName());
197     mgr->AddTask(taskPHOS);
198     if(fUseHLTTrigger==kFALSE)
199        AliAnalysisDataContainer *coutputPHOS =  mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison.root");  
200     else 
201        AliAnalysisDataContainer *coutputPHOS =  mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison_triggered.root");      
202     mgr->ConnectInput(taskPHOS,0,mgr->GetCommonInputContainer());
203     mgr->ConnectOutput(taskPHOS,1,coutputPHOS);
204   }
205
206   if(bEMCAL){
207     AliAnalysisTaskHLTEMCAL *taskEMCAL = new AliAnalysisTaskHLTEMCAL("offhlt_comparison_EMCAL");
208     taskEMCAL->SetUseHLTTriggerDecision(fUseHLTTrigger);
209     if(fUseHLTTrigger==kTRUE) printf("\n\nOnly HLT triggered events will be used to fill the distributions for task %s.\n\n", taskEMCAL->GetName());
210     mgr->AddTask(taskEMCAL);
211     if(fUseHLTTrigger==kFALSE)
212        AliAnalysisDataContainer *coutputEMCAL =  mgr->CreateContainer("emcal_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-EMCAL-comparison.root");  
213     else
214       AliAnalysisDataContainer *coutputEMCAL =  mgr->CreateContainer("emcal_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-EMCAL-comparison_triggered.root");      
215     mgr->ConnectInput(taskEMCAL,0,mgr->GetCommonInputContainer());
216     mgr->ConnectOutput(taskEMCAL,1,coutputEMCAL);
217   }
218   
219   if(bGLOBAL){
220     AliAnalysisTaskHLT *taskGLOBAL = new AliAnalysisTaskHLT("offhlt_comparison_GLOBAL");
221     taskGLOBAL->SetUseHLTTriggerDecision(fUseHLTTrigger);
222     if(fUseHLTTrigger==kTRUE) printf("\n\nOnly HLT triggered events will be used to fill the distributions for task %s.\n\n", taskGLOBAL->GetName());
223     //taskGLOBAL->SelectCollisionCandidates();
224     mgr->AddTask(taskGLOBAL);
225     if(fUseHLTTrigger==kFALSE)
226       AliAnalysisDataContainer *coutputGLOBAL =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison.root");  
227     else 
228       AliAnalysisDataContainer *coutputGLOBAL =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer,"HLT-OFFLINE-GLOBAL-comparison_triggered.root");  
229     mgr->ConnectInput(taskGLOBAL,0,mgr->GetCommonInputContainer());
230     mgr->ConnectOutput(taskGLOBAL,1,coutputGLOBAL);
231   }
232
233   if(bPWG1){
234     Bool_t hasMC=kFALSE;  
235     // -- Add Task for HLT and Offline
236     AliPerformanceTask *HLTtpcQA = AddTaskPerformance(hasMC,kFALSE,kTRUE);
237     AliPerformanceTask *tpcQA = AddTaskPerformance(hasMC,kFALSE); 
238     if(!HLTtpcQA || !tpcQA) {
239       Error("RunPerformanceTrain","AliPerformanceTask not created!");
240       return;
241     }
242   }
243   if(bD0){
244     float cuts[7]={0.5,0.04,0.7,0.8,0.05,-0.00025,0.7};
245     AliAnalysisTaskD0Trigger *taskD0 = new AliAnalysisTaskD0Trigger("offhlt_comparison_D0",cuts);
246     mgr->AddTask(taskD0);
247     AliAnalysisDataContainer *coutputD0 =  mgr->CreateContainer("D0_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-D0-comparison.root");  
248     mgr->ConnectInput(taskD0,0,mgr->GetCommonInputContainer());
249     mgr->ConnectOutput(taskD0,1,coutputD0);
250   }  
251   
252   if(bCB){
253      AliAnalysisTaskHLTCentralBarrel *taskCB = new AliAnalysisTaskHLTCentralBarrel("offhlt_comparison_CB"); 
254      mgr->AddTask(taskCB); 
255      taskCB->SetBeamType(beamType);
256      if(beamType.Contains("Pb-Pb")){
257         gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskCentrality.C");
258         AliCentralitySelectionTask *taskCentrality = AddTaskCentrality(); 
259         taskCentrality->SetPass(1);
260      }        
261      taskCB->SetOptions(options);
262      AliAnalysisDataContainer *coutputCB =  mgr->CreateContainer("esd_thnsparse", TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-CentralBarrel-comparison.root");  
263      mgr->ConnectInput(taskCB,0,mgr->GetCommonInputContainer());
264      mgr->ConnectOutput(taskCB,1,coutputCB);
265   }
266   
267   if (!mgr->InitAnalysis()) return;
268   mgr->PrintStatus();
269   mgr->StartAnalysis("local",chain, nEvents);
270
271   timer.Stop();
272   timer.Print();
273 }
274
275 void compare_HLT_offline_local(){
276   cout << " " << endl;
277   cout << " Usage examples:" << endl;
278   cout << "    compare-HLT-offline-local.C'(file, taskOption, taskFolder, beamType, options, fUseHLTTrigger, nEvents)' 2>&1 | tee log" << endl;
279   cout << "    compare-HLT-offline-local.C'(\"AliESDs.root\",\"global\")' 2>&1 | tee log" << endl;
280   cout << "    compare-HLT-offline-local.C'(\"AliESDs.root\",\"global\",\"./\", \"p-p\", \"event-off event-hlt track-off track-hlt\", kFALSE, nEvents)' 2>&1 | tee log" << endl;
281   cout << "    compare-HLT-offline-local.C'(\"AliESDs.root\",\"global phos cb D0\", \"./\", \"Pb-Pb\", \"event-hlt\", kTRUE, nEvents)' 2>&1 | tee log" << endl;
282   cout << "    compare-HLT-offline-local.C'(\"alien:///alice/data/2010/LHC10b/000115322/ESDs/pass1/10000115322040.20/AliESDs.root\",\"global\")' 2>&1 | tee log" << endl;
283   cout << " " << endl;
284 }