]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/tasks/macros/compare-HLT-offline-grid.C
- added argument to select the task location
[u/mrichter/AliRoot.git] / HLT / QA / tasks / macros / compare-HLT-offline-grid.C
1 // $Id$
2 /*
3  * Example macro to run an HLT analysis task for comparing the offline
4  * with the HLT esd tree on the GRID.
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 -q compare-HLT-offline-grid.C'("000115322","/alice/data/2010/LHC10b","ESDcomparison","output","full","global")'
12  * </pre>
13  * - run number
14  * - GRID input directory, where you define in which LHC period the run number belongs to
15  * - GRID working directory, where the .xml, .jdl and the task are uploaded (you have to create it yourself in advance)
16  * - GRID output directory with respect to the working one, where the output files of the task are located (you have to create it yourself in advance)
17  * - run in full mode, i.e. completely on the GRID with all the chunks of the run processed
18  * - specify the analysis task you want to run
19  * - specify whether you are interested only in HLT triggered events
20  *
21  * @ingroup alihlt_qa
22  * @author Hege.Erdal@student.uib.no, Kalliopi.Kanaki@ift.uib.no
23  */
24
25 void compare_HLT_offline_grid(TString runNumber, 
26                               TString dataDir, 
27                               TString gridWorkingDir, 
28                               TString gridOutputDir, 
29                               const char* mode = "full", 
30                               const char* detectorTask="global",
31                               TString taskFolder="$ALICE_ROOT/HLT/QA/tasks/",
32                               bool fUseHLTTrigger=kFALSE
33                              )
34 {
35  
36   TStopwatch timer;
37   timer.Start();
38
39   //gSystem->Load("libTree");
40   gSystem->Load("libCore");
41   gSystem->Load("libGeom");
42   gSystem->Load("libVMC");
43   gSystem->Load("libPhysics");
44  
45   //----------- Loading the required libraries ---------//
46
47   gSystem->Load("libSTEERBase");
48   gSystem->Load("libESD");
49   gSystem->Load("libAOD");
50   gSystem->Load("libANALYSIS");
51   gSystem->Load("libANALYSISalice");
52   gSystem->Load("libHLTbase");
53   gROOT->ProcessLine(".include $ALICE_ROOT/include"); 
54   gSystem->AddIncludePath("-I$ALICE_ROOT/HLT/BASE -I.")
55   //gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
56
57   
58   Bool_t bAll=kFALSE, bTPC=kFALSE, bPHOS=kFALSE, bEMCAL=kFALSE, bITS=kFALSE, bGLOBAL=kFALSE, bD0=kFALSE, bCB=kFALSE;
59  
60   TString allArgs = detectorTask;
61   TString argument;
62  
63   TObjArray *pTokens = allArgs.Tokenize(" ");
64   if(pTokens){
65      for(int i=0; i<pTokens->GetEntries(); i++){
66          argument=((TObjString*)pTokens->At(i))->GetString();
67          if(argument.IsNull()) continue;
68
69          if(argument.CompareTo("phos", TString::kIgnoreCase)==0){
70             bPHOS = kTRUE;
71             continue;
72          }
73          if(argument.CompareTo("emcal", TString::kIgnoreCase)==0){
74             bEMCAL = kTRUE;
75             continue;
76          }         
77          if(argument.CompareTo("global", TString::kIgnoreCase)==0){
78             bGLOBAL = kTRUE;
79             continue;
80          }  
81          if(argument.CompareTo("D0", TString::kIgnoreCase)==0){
82            bD0 = kTRUE;
83            continue;
84          }  
85          if(argument.CompareTo("cb", TString::kIgnoreCase)==0){
86            bCB = kTRUE;
87            continue;
88          }  
89          if(argument.CompareTo("all",TString::kIgnoreCase)==0){
90             bPHOS   = kTRUE;
91             bEMCAL  = kTRUE;
92             bGLOBAL = kTRUE;
93             bAll    = kTRUE;
94             continue;
95          }
96          else break;
97     }
98   }
99     
100   //-------- Make the analysis manager ---------------//
101  
102   AliAnalysisManager *mgr  = new AliAnalysisManager("TestManager");
103   AliESDInputHandler *esdH = new AliESDInputHandler;
104   esdH->SetReadHLT();
105   esdH->SetReadFriends(kFALSE);
106   mgr->SetInputEventHandler(esdH);  
107   mgr->SetNSysInfo(1000);
108
109   //To use Physics Selection
110   //AliPhysicsSelectionTask *physSelTask = AddTaskPhysicsSelection(kFALSE,kTRUE);
111
112   // Create and configure the alien handler plugin
113   gROOT->LoadMacro("CreateAlienHandler.C");
114   AliAnalysisGrid *alienHandler = CreateAlienHandler(runNumber, dataDir, gridWorkingDir, gridOutputDir, mode, detectorTask);
115   if (!alienHandler) return;
116   
117   // Connect plugin to the analysis manager
118   mgr->SetGridHandler(alienHandler);
119  
120   //-------------- Compile the analysis tasks ---------- //
121   if(bPHOS && bEMCAL) {
122     gSystem->Load("libHLTbase");
123     gSystem->Load("libAliHLTUtil");
124     gSystem->Load("libAliHLTGlobal");
125     TString strTask1("AliAnalysisTaskHLTCalo.cxx+");
126     TString strTask2("AliAnalysisTaskHLTPHOS.cxx+");    
127     TString strTask3("AliAnalysisTaskHLTEMCAL.cxx+");    
128     gROOT->LoadMacro(taskFolder+strTask1); 
129     gROOT->LoadMacro(taskFolder+strTask2);  
130     gROOT->LoadMacro(taskFolder+strTask3);      
131   }
132   else if(bPHOS) {
133     gSystem->Load("libHLTbase");
134     gSystem->Load("libAliHLTUtil");
135     gSystem->Load("libAliHLTGlobal");
136     TString strTask1("AliAnalysisTaskHLTCalo.cxx+");
137     TString strTask2("AliAnalysisTaskHLTPHOS.cxx+");    
138     gROOT->LoadMacro(taskFolder+strTask1); 
139     gROOT->LoadMacro(taskFolder+strTask2);  
140   }
141   else if(bEMCAL) {
142     gSystem->Load("libHLTbase");
143     gSystem->Load("libAliHLTUtil");
144     gSystem->Load("libAliHLTGlobal");
145     TString strTask1("AliAnalysisTaskHLTCalo.cxx+");
146     TString strTask2("AliAnalysisTaskHLTEMCAL.cxx+");    
147     gROOT->LoadMacro(taskFolder+strTask1); 
148     gROOT->LoadMacro(taskFolder+strTask2);  
149   }
150   if(bGLOBAL){
151      TString strTask("AliAnalysisTaskHLT.cxx+");   
152      gROOT->LoadMacro(taskFolder+strTask);
153   }
154   if(bD0){
155     TString strTask("AliAnalysisTaskD0Trigger.cxx+");
156     gROOT->LoadMacro(taskFolder+strTask);     
157   }
158   if(bCB){
159     TString strTask("AliAnalysisTaskHLTCentralBarrel.cxx+");
160     gROOT->LoadMacro(taskFolder+strTask);
161   } 
162  
163   //-------------- define the tasks ------------//
164   
165   if(bPHOS){
166      AliAnalysisTaskHLTPHOS *taskPHOS = new AliAnalysisTaskHLTPHOS("offhlt_comparison_PHOS");
167      mgr->AddTask(taskPHOS);
168      AliAnalysisDataContainer *coutputPHOS =  mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison.root");  
169      mgr->ConnectInput(taskPHOS,0,mgr->GetCommonInputContainer());
170      mgr->ConnectOutput(taskPHOS,1,coutputPHOS);
171   }
172   if(bEMCAL){
173      AliAnalysisTaskHLTEMCAL *taskEMCAL = new AliAnalysisTaskHLTEMCAL("offhlt_comparison_EMCAL");
174      mgr->AddTask(taskEMCAL);
175      AliAnalysisDataContainer *coutputEMCAL =  mgr->CreateContainer("emcal_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-EMCAL-comparison.root");  
176      mgr->ConnectInput(taskEMCAL,0,mgr->GetCommonInputContainer());
177      mgr->ConnectOutput(taskEMCAL,1,coutputEMCAL);
178   }
179   
180   if(bGLOBAL){
181      AliAnalysisTaskHLT *taskGLOBAL = new AliAnalysisTaskHLT("offhlt_comparison_GLOBAL");
182      taskGLOBAL->SetUseHLTTriggerDecision(fUseHLTTrigger);
183      if(fUseHLTTrigger==kTRUE) printf("\n\nOnly HLT triggered events will be used to fill the distributions for task %s.\n\n", taskGLOBAL->GetName());
184      //taskGLOBAL->SelectCollisionCandidates();
185      mgr->AddTask(taskGLOBAL);
186      if(fUseHLTTrigger==kFALSE)
187        AliAnalysisDataContainer *coutputGLOBAL =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison.root");  
188      else
189        AliAnalysisDataContainer *coutputGLOBAL =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison_triggered.root");  
190      mgr->ConnectInput(taskGLOBAL,0,mgr->GetCommonInputContainer());
191      mgr->ConnectOutput(taskGLOBAL,1,coutputGLOBAL);
192   }
193   if(bD0){
194     float cuts[7]={0.5,0.04,0.7,0.8,0.05,-0.00025,0.7};
195     AliAnalysisTaskD0Trigger *taskD0 = new AliAnalysisTaskD0Trigger("offhlt_comparison_D0_Trigger",cuts);
196     mgr->AddTask(taskD0);
197     AliAnalysisDataContainer *coutputD0 =  mgr->CreateContainer("D0_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-D0-comparison.root");  
198     mgr->ConnectInput(taskD0,0,mgr->GetCommonInputContainer());
199     mgr->ConnectOutput(taskD0,1,coutputD0);
200   }
201   if(bCB){
202      AliAnalysisTaskHLTCentralBarrel *taskCB = new AliAnalysisTaskHLTCentralBarrel("offhlt_comparison_CB");
203      mgr->AddTask(taskCB);
204      AliAnalysisDataContainer *coutputCB =  mgr->CreateContainer("esd_thnsparse",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-CentralBarrel-comparison.root");       
205      mgr->ConnectInput(taskCB,0,mgr->GetCommonInputContainer());
206      mgr->ConnectOutput(taskCB,1,coutputCB);
207   }
208   // Enable debug printouts
209   mgr->SetDebugLevel(2);
210
211   if (!mgr->InitAnalysis()) return;
212   mgr->PrintStatus();
213   mgr->StartAnalysis("grid");
214
215   timer.Stop();
216   timer.Print();
217 }
218
219 void compare_HLT_offline_grid(){
220   cout << " " << endl;
221   cout << " Usage examples:" << endl;
222   cout << "    compare-HLT-offline-grid.C'(runNumber, dataDir, gridWorkingDir, gridOutputDir, mode, taskOption, taskFolder, fUseHLTTrigger)' 2>&1 | tee log" << endl;
223   cout << "    compare-HLT-offline-grid.C'(\"000115322\",\"/alice/data/2010/LHC10b\",\"ESDcomparison\",\"output\",\"full\",\"global\",\"./\",kTRUE)' 2>&1 | tee log" << endl;
224   cout << " " << endl;
225 }