]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/QA/tasks/macros/compare-HLT-offline-local.C
- added check for the existence of HLT tracks in the HLTesdTree and printout of the...
[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'("/home/blabla/AliESDs.root","phos")' 2>&1 | tee task.log
12  *   aliroot -b -q -l compare_HLT_offline_local.C'("/home/blabla/AliESDs.root","phos tpc global")' 2>&1 | tee task.log
13  *   aliroot -q compare-HLT-offline-local.C'("alien:///alice/data/2010/LHC10b/000115322/ESDs/pass1/10000115322040.20/AliESDs.root","global")' 2>&1 | tee log
14  * </pre>
15  * 
16  * If alien:// is contained in the name of the file, then the macro connects to the grid to access the file.
17  *
18  * @ingroup alihlt_qa
19  * @author zbyin@mail.ccnu.edu.cn, Kalliopi.Kanaki@ift.uib.no
20  */
21
22 void compare_HLT_offline_local(TString file, 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   if(!AliAnalysisGrid::CreateToken()) return NULL;
89   
90   if(file.Contains("alien")) TGrid::Connect("alien://");
91
92   TChain *chain = new TChain("esdTree"); 
93   chain->Add(file);
94  
95   //chain->Add("alien:///alice/data/2010/LHC10b/000115322/ESDs/pass1/10000115322040.20/AliESDs.root");
96   //chain->Add("...");
97    
98   //-------- Make the analysis manager ---------------//
99  
100   AliAnalysisManager *mgr  = new AliAnalysisManager("TestManager");
101   AliESDInputHandler *esdH = new AliESDInputHandler;
102   esdH->SetReadHLT();
103   esdH->SetReadFriends(kFALSE);
104   mgr->SetInputEventHandler(esdH);  
105   mgr->SetNSysInfo(1000);
106  
107   //-------------- define the tasks ------------//
108   
109   if(bTPC){ 
110      AliAnalysisTaskHLTTPC *taskTPC = new AliAnalysisTaskHLTTPC("offhlt_comparison_TPC");
111      mgr->AddTask(taskTPC);
112      AliAnalysisDataContainer *coutput1 =  mgr->CreateContainer("tpc_histograms", TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-TPC-comparison.root");  
113      mgr->ConnectInput(taskTPC,0,mgr->GetCommonInputContainer());
114      mgr->ConnectOutput(taskTPC,1,coutput1);
115   }
116
117   if(bPHOS){
118      AliAnalysisTaskHLTPHOS *taskPHOS = new AliAnalysisTaskHLTPHOS("offhlt_comparison_PHOS");
119      mgr->AddTask(taskPHOS);
120      AliAnalysisDataContainer *coutput2 =  mgr->CreateContainer("phos_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-PHOS-comparison.root");  
121      mgr->ConnectInput(taskPHOS,0,mgr->GetCommonInputContainer());
122      mgr->ConnectOutput(taskPHOS,1,coutput2);
123   }
124   
125   if(bITS){
126      AliAnalysisTaskHLTITS *taskITS = new AliAnalysisTaskHLTITS("offhlt_comparison_ITS");
127      mgr->AddTask(taskITS);
128      AliAnalysisDataContainer *coutput3 =  mgr->CreateContainer("its_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-ITS-comparison.root");  
129      mgr->ConnectInput(taskITS,0,mgr->GetCommonInputContainer());
130      mgr->ConnectOutput(taskITS,1,coutput3);
131   }
132  
133   if(bGLOBAL){
134      AliAnalysisTaskHLT *taskGLOBAL = new AliAnalysisTaskHLT("offhlt_comparison_GLOBAL");
135      mgr->AddTask(taskGLOBAL);
136      AliAnalysisDataContainer *coutput4 =  mgr->CreateContainer("global_histograms",TList::Class(), AliAnalysisManager::kOutputContainer, "HLT-OFFLINE-GLOBAL-comparison.root");  
137      mgr->ConnectInput(taskGLOBAL,0,mgr->GetCommonInputContainer());
138      mgr->ConnectOutput(taskGLOBAL,1,coutput4);
139   }
140   
141   if (!mgr->InitAnalysis()) return;
142   mgr->PrintStatus();
143   mgr->StartAnalysis("local",chain);
144
145   timer.Stop();
146   timer.Print();
147 }