]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/macros/lego_train/RunALICE.C
PWG2/SPECTRA -> PWGLF/SPECTRA migration
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / macros / lego_train / RunALICE.C
1 #ifndef __CINT__
2 #include <TSystem.h>
3 #include <TROOT.h>
4 #include <Rtypes.h>
5 #include <TString.h>
6 #include <TNamed.h>
7 #include <TObjArray.h>
8 #include <TObjString.h>
9 #include <TList.h>
10 #include <ANALYSIS/AliAnalysisManager.h>
11 #include <ANALYSIS/AliMultiInputEventHandler.h>
12 #include <STEER/AOD/AliAODHandler.h>
13 #include <TStopwatch.h>
14 #endif
15
16 Bool_t RunALICE(TString anSrc = "grid",
17                 TString anMode = "terminate",
18                 TString input="aod" /*or "esd"*/,
19                 TString inputMC="" /*or "mc"*/,
20                 Long64_t nEvents = 1e10,
21                 Long64_t nSkip = 0,
22                 TString dsName="",
23                 TString alirsnliteManagers ="",
24                 Bool_t useMultiHandler=kTRUE,
25                 TString alirsnlitesrc ="$ALICE_ROOT",
26                 TString alirsnlitetasks =""
27                ) {
28
29    // some init work
30    anSrc.ToLower(); anMode.ToLower(); input.ToLower(); inputMC.ToLower();
31
32    // loads libs and setup include paths
33    if (LoadLibsBase(alirsnlitesrc)) return kFALSE;
34
35    // reset manager if already exists
36    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
37    if (mgr) delete mgr;
38    mgr = new AliAnalysisManager("AliRsnLiteAM","AliRsnLite Analysis Manager");
39
40    Bool_t useAODOut = kFALSE;
41    CreateInputHandlers(input,inputMC,useAODOut,useMultiHandler);
42
43    // add default grid handler
44    gROOT->LoadMacro("SetupAnalysisPlugin.C");
45    AliAnalysisGrid *analysisPlugin = SetupAnalysisPlugin(anMode.Data());
46    if (!analysisPlugin) { Printf("Error : analysisPlugin is null !!!"); return kFALSE; }
47    mgr->SetGridHandler(analysisPlugin);
48    if (!dsName.IsNull()) {
49       if (!anSrc.CompareTo("proof") && !anMode.CompareTo("full")) {
50          analysisPlugin->SetProofDataSet(dsName.Data());
51          Printf(Form("Using DataSet %s ...",dsName.Data()));
52       } else {
53          analysisPlugin->SetFileForTestMode(dsName.Data());
54          Printf(Form("Using Test file %s ...",dsName.Data()));
55       }
56    }
57
58    TList *listManagers = CreateListOfManagersFromDir(alirsnliteManagers,alirsnlitetasks);
59    if (!listManagers) { Printf("Error : CreateListOfManagersFromDir failed !!!"); return kFALSE;}
60
61    // adds all tasks
62    if (!AddAllManagers(listManagers, anSrc, anMode,input,inputMC)) { Printf("Error : AddAllManagers failed !!!"); return kFALSE;}
63
64    TStopwatch timer;
65    timer.Start();
66    // runs analysis
67    if (!RunAnalysisManager(anSrc, anMode.Data(), nEvents, nSkip)) { Printf("Error : RunAnalysisManager failed !!!"); return kFALSE;}
68
69    timer.Stop();
70    timer.Print();
71    Printf("Working directory is %s ...", gSystem->WorkingDirectory());
72    Printf("Done OK");
73    return kTRUE;
74
75 }
76
77 Int_t LoadLibsBase(TString alirsnlitesrc) {
78    Int_t num = 0;
79    if (gSystem->Load("libTree.so") < 0) {num++; return num;}
80    if (gSystem->Load("libGeom.so") < 0) {num++; return num;}
81    if (gSystem->Load("libVMC.so") < 0) {num++; return num;}
82    if (gSystem->Load("libMinuit.so") < 0) {num++; return num;}
83    if (gSystem->Load("libPhysics.so") < 0) {num++; return num;}
84    if (gSystem->Load("libSTEERBase.so") < 0) {num++; return num;}
85    if (gSystem->Load("libESD.so") < 0) {num++; return num;}
86    if (gSystem->Load("libAOD.so") < 0) {num++; return num;}
87    if (gSystem->Load("libANALYSIS.so") < 0) {num++; return num;}
88    if (gSystem->Load("libOADB.so") < 0) {num++; return num;}
89    if (gSystem->Load("libANALYSISalice.so") < 0) {num++; return num;}
90
91    gSystem->AddIncludePath(Form("-I\"%s/include\"", gSystem->ExpandPathName(alirsnlitesrc.Data())));
92    gROOT->ProcessLine(Form(".include %s/include", gSystem->ExpandPathName(alirsnlitesrc.Data())));
93 }
94
95 Bool_t CreateInputHandlers(TString input,TString inputMC,Bool_t useAODOut=kFALSE,Bool_t useMultiHandler=kTRUE) {
96
97    Bool_t useMC = !inputMC.CompareTo("mc");
98
99    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
100    if (!mgr) { Printf("Error [CreateInputHandlers] : mgr is null !!!"); return kFALSE; }
101
102    if (useMultiHandler) {
103       AliMultiInputEventHandler *inputHandler = new AliMultiInputEventHandler();
104       if (!input.CompareTo("esd")) {
105          inputHandler->AddInputEventHandler(new AliESDInputHandler());
106          if (useMC) inputHandler->AddInputEventHandler(new AliMCEventHandler());
107       } else if (!input.CompareTo("aod")) {
108          inputHandler->AddInputEventHandler(new AliAODInputHandler());
109       }
110
111       mgr->SetInputEventHandler(inputHandler);
112    } else {
113       if (!input.CompareTo("esd")) {
114          mgr->SetInputEventHandler(new AliESDInputHandler());
115          if (useMC) mgr->SetMCtruthEventHandler(new AliMCEventHandler());
116       } else if (!input.CompareTo("aod")) {
117          mgr->SetInputEventHandler(new AliAODInputHandler());
118       }
119       mgr->SetInputEventHandler(inputHandler);
120    }
121
122    if (useAODOut) {
123       AliAODHandler *aodHandler   = new AliAODHandler();
124       aodHandler->SetOutputFileName("AliAOD.root");
125       mgr->SetOutputEventHandler(aodHandler);
126    }
127
128 }
129
130 TList *CreateListOfManagersFromDir(TString listManagersNames="",TString dir="") {
131
132    TList *listManagers = new TList;
133    TString dirsStr;
134    TObjArray *dirs=0;
135
136    if (listManagersNames.IsNull()) {
137       if (dir.IsNull() || gSystem->AccessPathName(gSystem->ExpandPathName(dir.Data()))) {
138          Printf("Error [CreateListOfManagersFromDir] : Dir '%s' doesn't exists !!!",dir.Data());
139          return 0;
140       }
141       dirsStr = gSystem->GetFromPipe(Form("ls %s",dir.Data()));
142       dirs = dirsStr.Tokenize("\n");
143    } else {
144       dirsStr = listManagersNames;
145       dirs = dirsStr.Tokenize(" ");
146    }
147
148    TIter next(dirs);
149    Int_t counter=0;
150    TObjString *str,*strtmp;
151    TObjArray *mydirstrTok;
152    TString mydirstr,main,prefix;
153    while ((str = (TObjString *)next.Next())) {
154       // TODO add direcotry
155       mydirstr = str->GetString();
156       if (mydirstr.IsNull()) continue;
157
158       Printf("Adding %s .,,",mydirstr.Data());
159       mydirstrTok = mydirstr.Tokenize("_");
160
161       main = ((TObjString *)mydirstrTok->At(0))->GetString();
162
163       strtmp = (TObjString *)mydirstrTok->At(1);
164       if (strtmp) prefix = strtmp->GetString(); else prefix="";
165
166       listManagers->Add(new TNamed(main,prefix));
167    }
168
169    return listManagers;
170 }
171
172 Bool_t AddAllManagers(TList *listManagers,TString anSrc, TString anMode,TString input,TString inputMC) {
173    TIter next(listManagers);
174    Int_t counter=0;
175    TNamed *name;
176    while ((name = (TNamed *)next.Next())) {
177       if (!AddAnalysisManager(name->GetName(), anSrc, anMode,input,inputMC,name->GetTitle(),Form("%d",counter++))) {
178          Printf("Error: Problem adding %s",name->GetName());
179          return kFALSE;
180       }
181    }
182
183    return kTRUE;
184 }
185
186 Bool_t AddAnalysisManager(TString managerMacro, TString anSrc, TString anMode,TString input,TString inputMC, TString postfix,TString idStr) {
187    gROOT->LoadMacro(Form("%s.C", managerMacro.Data()));
188    return gROOT->ProcessLine(Form("%s\(\"%s\",\"%s\",\"%s\"\,\"%s\",\"%s\",\"%s\"\);", managerMacro.Data(), anSrc.Data(), anMode.Data(),input.Data(),inputMC.Data(), postfix.Data(),idStr.Data()));
189 }
190
191 Bool_t RunAnalysisManager(TString anSrc = "proof", TString anMode = "test", Long64_t nEvents = 1e10, Long64_t nSkip = 0) {
192
193    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
194
195    if (!mgr) { Printf("Error [RunAnalysisManager] : mgr is null !!!"); return kFALSE; }
196
197    // Run analysis
198    mgr->InitAnalysis();
199    mgr->PrintStatus();
200
201    if ((!anSrc.CompareTo("proof")) || (!anSrc.CompareTo("local"))) {
202       mgr->StartAnalysis(anSrc.Data(), nEvents, nSkip);
203    } else {
204       mgr->StartAnalysis(anSrc.Data());
205    }
206
207    return kTRUE;
208 }