]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/tasks/macros/compare-HLT-offline-local.C
- implemented function that allows selection of HLT triggered events (Hege)
[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  * Usage:
12  * <pre>
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
16  * </pre>
17  * 
18  * If alien:// is contained in the name of the file, then the macro connects to the grid to access the file.
19  * 
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.
22  *
23  * @ingroup alihlt_qa
24  * @author Kalliopi.Kanaki@ift.uib.no, Hege.Erdal@student.uib.no
25  */
26
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
32                               )
33 {
34
35   TStopwatch timer;
36   timer.Start();
37
38   gSystem->Load("libTree.so");
39   gSystem->Load("libGeom.so");
40   gSystem->Load("libVMC.so");
41   gSystem->Load("libPhysics.so");
42  
43   //----------- Loading the required libraries ---------//
44
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");
51  
52   gSystem->AddIncludePath("-I$ALICE_ROOT/PWG1/TPC -I.");
53   
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");
61  
62   gROOT->ProcessLine(".include $ALICE_ROOT/include");
63   //gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
64   
65   Bool_t bPHOS=kFALSE, bGLOBAL=kFALSE, bEMCAL = kFALSE, bPWG1 = kFALSE, bD0=kFALSE;
66  
67   TString allArgs = detectorTask;
68   TString argument;
69  
70   TObjArray *pTokens = allArgs.Tokenize(" ");
71   if(pTokens){
72     for(int i=0; i<pTokens->GetEntries(); i++){
73       argument=((TObjString*)pTokens->At(i))->GetString();
74       if(argument.IsNull()) continue;
75
76       if(argument.CompareTo("phos", TString::kIgnoreCase)==0){
77         bPHOS = kTRUE;
78         continue;
79       }         
80       else if(argument.CompareTo("emcal", TString::kIgnoreCase)==0){
81         bEMCAL = kTRUE;
82         continue;
83       }         
84       if(argument.CompareTo("global", TString::kIgnoreCase)==0){
85         bGLOBAL = kTRUE;
86         continue;
87       }      
88       if(argument.CompareTo("pwg1", TString::kIgnoreCase)==0){
89         bPWG1 = kTRUE;
90         continue;
91       }
92       if(argument.CompareTo("D0", TString::kIgnoreCase)==0){
93         bD0 = kTRUE;
94         continue;
95       }
96       if(argument.CompareTo("all",TString::kIgnoreCase)==0){
97         bPHOS   = kTRUE;
98         bEMCAL  = kTRUE;
99         bGLOBAL = kTRUE; 
100         bD0     = kTRUE;   
101         continue;
102       }
103       else break;
104     }
105   }
106       
107   //-------------- Compile the analysis tasks ---------- //
108   
109   if(bPHOS){
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;
118   }
119   
120   if(bEMCAL){
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;
129   }  
130   
131   if(bGLOBAL){
132      TString strTask("AliAnalysisTaskHLT.cxx+");
133      gROOT->LoadMacro(taskFolder+strTask);
134      cout << "\n========= You are loading the following task --> "<< taskFolder+strTask  << endl;
135   }
136   if(bD0){
137      TString strTask("AliAnalysisTaskD0Trigger.cxx+");
138      gROOT->LoadMacro(taskFolder+strTask); 
139      cout << "\n========= You are loading the following task --> "<< taskFolder+strTask  << endl;
140   }
141   
142   if(bPWG1) gROOT->LoadMacro("$ALICE_ROOT/HLT/QA/tasks/macros/AddTaskPerformance.C");
143    
144   if(file.Contains("alien")) TGrid::Connect("alien://");
145     
146   if(file.Contains("AliESDs.root")){
147     TChain *chain = new TChain("esdTree"); 
148     chain->Add(file);
149   }
150   
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());
156   }
157
158   if(!chain){
159     Printf("Chain is empty");
160     return;
161   }
162
163   //To only select HLT triggered events
164   //Bool_t fUseHLTTrigger=kFALSE;
165    
166   //-------- Make the analysis manager ---------------//
167  
168   AliAnalysisManager *mgr  = new AliAnalysisManager("TestManager");
169   AliESDInputHandler *esdH = new AliESDInputHandler;
170
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);
176
177   //AliPhysicsSelectionTask *physSelTask = AddTaskPhysicsSelection(kFALSE,kTRUE);
178  
179   //-------------- define the tasks ------------//
180   
181   if(bPHOS){
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);
187   }
188
189   if(bEMCAL){
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);
195   }
196   
197   if(bGLOBAL){
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");  
205     else 
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);
209   }
210
211   if(bPWG1){
212     Bool_t hasMC=kFALSE;  
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!");
218       return;
219     }
220   }
221   if(bD0){
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);
228   }
229   
230   if (!mgr->InitAnalysis()) return;
231   mgr->PrintStatus();
232   mgr->StartAnalysis("local",chain, nEvents);
233
234   timer.Stop();
235   timer.Print();
236 }
237
238 void compare_HLT_offline_local(){
239   cout << " " << endl;
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;
247   cout << " " << endl;
248 }