]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/run.C
update for the new Rezidual plots
[u/mrichter/AliRoot.git] / TRD / qaRec / run.C
1 // Steer TRD QA train for Reconstruction (Clusterizer, Tracking and PID).
2 // 
3 // Usage:
4 //   run.C(tasks, files, entries)
5 //   tasks : "ALL" or one/more of the following:
6 //     "EFF"  : TRD Tracking Efficiency 
7 //     "EFFC" : TRD Tracking Efficiency Combined (barrel + stand alone) - only in case of simulations
8 //     "RES"  : TRD tracking Resolution
9 //     "CLRES": clusters Resolution
10 //     "CAL"  : TRD calibration
11 //     "PID"  : TRD PID - pion efficiency 
12 //     "PIDR" : TRD PID - reference data
13 //     "DET"  : Basic TRD Detector checks
14 //     "NOFR" : Data set does not have AliESDfriends.root 
15 //     "NOMC" : Data set does not have Monte Carlo Informations (real data), so all tasks which rely
16 //              on MC information are switched off
17 //
18 // In compiled mode : 
19 // Don't forget to load first the libraries
20 // gSystem->Load("libMemStat.so")
21 // gSystem->Load("libMemStatGui.so")
22 // gSystem->Load("libANALYSIS.so")
23 // gSystem->Load("libTRDqaRec.so")
24 //
25 // Authors:
26 //   Alex Bercuci (A.Bercuci@gsi.de) 
27 //   Markus Fasel (m.Fasel@gsi.de) 
28
29 #ifndef __CINT__
30 #include <Riostream.h>
31
32 #include "TStopwatch.h"
33 #include "TMemStat.h"
34 #include "TMemStatViewerGUI.h"
35
36 #include "TROOT.h"
37 #include "TSystem.h"
38 #include "TError.h"
39 #include "TChain.h"
40
41 #include "AliMagFMaps.h"
42 #include "AliTracker.h"
43 #include "AliLog.h"
44 #include "AliCDBManager.h"
45 #include "AliAnalysisManager.h"
46 #include "AliAnalysisDataContainer.h"
47 #include "AliMCEventHandler.h"
48 #include "AliESDInputHandler.h"
49
50 #include "TRD/AliTRDtrackerV1.h"
51 #include "TRD/AliTRDcalibDB.h"
52 #include "TRD/qaRec/AliTRDtrackInfo/AliTRDeventInfo.h"
53 #include "TRD/qaRec/AliTRDtrackInfoGen.h"
54 #include "TRD/qaRec/AliTRDtrackingEfficiency.h"
55 #include "TRD/qaRec/AliTRDtrackingEfficiencyCombined.h"
56 #include "TRD/qaRec/AliTRDtrackingResolution.h"
57 #include "TRD/qaRec/AliTRDcalibration.h"
58 #include "TRD/qaRec/AliTRDpidChecker.h"
59 #include "TRD/qaRec/AliTRDpidRefMaker.h"
60 #include "TRD/qaRec/AliTRDcheckDetector.h"
61 #include "TRD/qaRec/AliTRDclusterResolution.h"
62 #endif
63
64 #include "run.h"
65
66 Bool_t MEM = kFALSE;
67
68 TChain* CreateESDChain(const char* filename = 0x0, Int_t nfiles=-1 );
69 void run(Char_t *tasks="ALL", const Char_t *files=0x0, Int_t nmax=-1)
70 {
71   TMemStat *mem = 0x0;
72   if(MEM){ 
73     gSystem->Load("libMemStat.so");
74     gSystem->Load("libMemStatGui.so");
75     mem = new TMemStat("new, gnubuildin");
76     mem->AddStamp("Start");
77   }
78
79   TStopwatch timer;
80   timer.Start();
81
82   if(gSystem->Load("libANALYSIS.so")<0) return;
83   if(gSystem->Load("libTRDqaRec.so")<0) return;
84   
85   Bool_t fHasMCdata = kTRUE;
86   Bool_t fHasFriends = kTRUE;
87   TObjArray *tasksArray = TString(tasks).Tokenize(" ");
88
89   Int_t fSteerTask = 0; SETBIT(fSteerTask, kInfoGen);
90   for(Int_t isel = 0; isel < tasksArray->GetEntriesFast(); isel++){
91     TString s = (dynamic_cast<TObjString *>(tasksArray->UncheckedAt(isel)))->String();
92     if(s.CompareTo("ALL") == 0){
93       for(Int_t itask = 1; itask <= NTRDTASKS; itask++) SETBIT(fSteerTask, itask);
94       continue;
95     } else if(s.CompareTo("NOFR") == 0){ 
96       fHasFriends = kFALSE;
97     } else if(s.CompareTo("NOMC") == 0){ 
98       fHasMCdata = kFALSE;
99     } else { 
100       Bool_t foundOpt = kFALSE;  
101       for(Int_t itask = 1; itask <= NTRDTASKS; itask++){
102         if(s.CompareTo(fgkTRDtaskOpt[itask]) != 0) continue;
103         SETBIT(fSteerTask, itask);
104         foundOpt = kTRUE;
105         break;
106       }
107       if(!foundOpt) Info("run.C", Form("Task %s not implemented (yet).", s.Data()));
108     }
109   }
110   // extra rules for calibration tasks
111   if(TSTBIT(fSteerTask, kClusterErrorParam)) SETBIT(fSteerTask, kTrackingResolution);
112
113   // define task list pointers;
114   AliTRDrecoTask *taskPtr[NTRDTASKS], *task = 0x0;
115   memset(taskPtr, 0, NTRDTASKS*sizeof(AliAnalysisTask*));
116
117   //____________________________________________//
118   //gROOT->LoadMacro(Form("%s/TRD/qaRec/CreateESDChain.C", gSystem->ExpandPathName("$ALICE_ROOT")));
119   TChain *chain = CreateESDChain(files, nmax);
120   //chain->SetBranchStatus("*", 0);
121   chain->SetBranchStatus("*FMD*",0);
122   chain->SetBranchStatus("*Calo*",0);
123   chain->SetBranchStatus("Tracks", 1);
124   chain->SetBranchStatus("ESDfriend*",1);
125   chain->Lookup();
126   chain->GetListOfFiles()->Print();
127   printf("\n ----> CHAIN HAS %d ENTRIES <----\n\n", (Int_t)chain->GetEntries());
128   
129   AliLog::SetGlobalLogLevel(AliLog::kError);
130
131   //____________________________________________
132   // Make the analysis manager
133   AliAnalysisManager *mgr = new AliAnalysisManager("TRD QA Reconstruction Manager");
134   //mgr->SetSpecialOutputLocation(source); // To Be Changed
135   AliVEventHandler *esdH = 0x0, *mcH = 0x0;
136   mgr->SetInputEventHandler(esdH = new AliESDInputHandler);
137   if(fHasMCdata) mgr->SetMCtruthEventHandler(mcH = new AliMCEventHandler());
138   //mgr->SetDebugLevel(10);
139
140   //____________________________________________
141   // TRD track summary generator
142   mgr->AddTask(task = new AliTRDtrackInfoGen());
143   taskPtr[(Int_t)kInfoGen] = task;
144   task->SetDebugLevel(1);
145   task->SetMCdata(fHasMCdata);
146   // Create containers for input/output
147   AliAnalysisDataContainer *cinput1 = mgr->CreateContainer("data", TChain::Class(), AliAnalysisManager::kInputContainer);
148   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("trackInfo", TObjArray::Class(), AliAnalysisManager::kExchangeContainer);
149   AliAnalysisDataContainer *coutput1a = mgr->CreateContainer("eventInfo", AliTRDeventInfo::Class(), AliAnalysisManager::kExchangeContainer);
150   mgr->ConnectInput( task, 0, cinput1);
151   mgr->ConnectOutput(task, 0, coutput1);
152   mgr->ConnectOutput(task, 1, coutput1a);
153
154   //____________________________________________
155   // TRD detector checker
156         if(TSTBIT(fSteerTask, kCheckDetector)){
157     mgr->AddTask(task = new AliTRDcheckDetector());
158     taskPtr[(Int_t)kCheckDetector] = task;
159     task->SetDebugLevel(4);
160     task->SetMCdata(fHasMCdata);
161     
162     // Create containers for input/output
163     mgr->ConnectInput( task, 0, coutput1);
164     mgr->ConnectInput( task, 1, coutput1a);
165     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
166   }
167
168   //____________________________________________
169   // TRD barrel tracking efficiency
170   if(fHasMCdata && TSTBIT(fSteerTask, kTrackingEfficiency)){
171     mgr->AddTask(task = new AliTRDtrackingEfficiency());
172     taskPtr[(Int_t)kTrackingEfficiency] = task;
173     task->SetDebugLevel(0);
174
175     //Create containers for input/output
176     mgr->ConnectInput( task, 0, coutput1);
177     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
178   }
179
180   //____________________________________________
181   // TRD combined tracking efficiency
182   if(fHasMCdata && TSTBIT(fSteerTask, kTrackingCombinedEfficiency)){
183     mgr->AddTask(task = new AliTRDtrackingEfficiencyCombined());
184     taskPtr[(Int_t)kTrackingCombinedEfficiency] = task;
185     task->SetDebugLevel(0);
186
187     // Create containers for input/output
188     mgr->ConnectInput( task, 0, coutput1);
189     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
190   }
191
192   //____________________________________________
193   // TRD tracking resolution
194   if(TSTBIT(fSteerTask, kTrackingResolution)){
195     mgr->AddTask(task = new AliTRDtrackingResolution());
196     taskPtr[(Int_t)kTrackingResolution] = task;
197     task->SetMCdata(fHasMCdata);
198     task->SetPostProcess(kFALSE);
199     task->SetDebugLevel(1);
200     
201     // Create containers for input/output
202     mgr->ConnectInput( task, 0, coutput1);
203     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
204
205     // Create output containers for calibration tasks
206     const Int_t nc = 5;
207     const Char_t *cn[nc] = {"ClRez", "TrkltRez", "TrkltPhiRez", "ClRes", "TrkltRes"}; 
208     AliAnalysisDataContainer *co[nc]; 
209     for(Int_t ic = 0; ic<nc; ic++){
210       co[ic] = mgr->CreateContainer(Form("%s%s", task->GetName(), cn[ic]), TObjArray::Class(), AliAnalysisManager::kExchangeContainer);
211       mgr->ConnectOutput(task, 1+ic, co[ic]);
212     }
213     
214     // test reconstruction calibration plugin
215     if(TSTBIT(fSteerTask, kClusterErrorParam)){
216       mgr->AddTask(task = new AliTRDclusterResolution());
217       taskPtr[(Int_t)kClusterErrorParam] = task;
218       mgr->ConnectInput(task, 0, co[0]);
219       mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
220   
221       mgr->AddTask(task = new AliTRDclusterResolution());
222       mgr->ConnectInput(task, 0, co[1]);
223       mgr->ConnectOutput(task, 0, mgr->CreateContainer(Form("%sMC", task->GetName()), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sMC.root", task->GetName())));
224     }
225   }
226
227   //____________________________________________
228   // TRD calibration
229   if(TSTBIT(fSteerTask, kCalibration)){
230     mgr->AddTask(task = new AliTRDcalibration());
231     taskPtr[(Int_t)kCalibration] = task;
232     ((AliTRDcalibration*)task)->SetLow(0);
233     ((AliTRDcalibration*)task)->SetHigh(30);
234     ((AliTRDcalibration*)task)->SetFillZero(kFALSE);
235     task->SetDebugLevel(0);
236
237     // Create containers for input/output
238     mgr->ConnectInput(task, 0, coutput1);
239     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
240   }
241   
242   //____________________________________________
243   // TRD pid checker
244   if(TSTBIT(fSteerTask, kPIDChecker)){
245     mgr->AddTask(task = new AliTRDpidChecker());
246     taskPtr[(Int_t)kPIDChecker] = task;
247     task->SetDebugLevel(0);
248     task->SetMCdata(fHasMCdata);
249     
250     // Create containers for input/output
251     mgr->ConnectInput( task, 0, coutput1);
252     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
253   }
254
255
256   //____________________________________________
257   // TRD pid reference 
258   if(TSTBIT(fSteerTask, kPIDRefMaker)){
259     mgr->AddTask(task = new AliTRDpidRefMaker());
260     taskPtr[(Int_t)kPIDRefMaker] = task;
261     task->SetDebugLevel(0);
262     task->SetMCdata(fHasMCdata);
263     
264     // Create containers for input/output
265     mgr->ConnectInput( task, 0, coutput1);
266     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
267     mgr->ConnectOutput(task, 1, mgr->CreateContainer(Form("%sNN", task->GetName()), TTree::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sNN.root", task->GetName())));
268     mgr->ConnectOutput(task, 2, mgr->CreateContainer(Form("%sLQ", task->GetName()), TTree::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sLQ.root", task->GetName())));
269   }
270
271
272   if (!mgr->InitAnalysis()) return;
273   printf("\n\tRUNNING TRAIN FOR TASKS:\n");
274   for(Int_t itask = 1; itask <= NTRDTASKS; itask++){
275     if(TSTBIT(fSteerTask, itask)) printf("\t   %s [%s]\n", taskPtr[itask]->GetTitle(), taskPtr[itask]->GetName());
276   }
277   printf("\n\n");
278   //mgr->PrintStatus();
279
280
281   AliCDBManager *cdbManager = AliCDBManager::Instance();
282   cdbManager->SetDefaultStorage("local://$ALICE_ROOT");
283   //cdbManager->SetSpecificStorage("TRD/Calib/FEE","local:///u/bailhach/aliroot/database30head/database");
284   cdbManager->SetRun(0);
285   cdbManager->SetCacheFlag(kFALSE);
286  
287   // initialize TRD settings
288   AliMagFMaps *field = new AliMagFMaps("Maps","Maps", 2, 1., 10., AliMagFMaps::k5kG);
289   AliTracker::SetFieldMap(field, kTRUE);
290   AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
291   AliTRDtrackerV1::SetNTimeBins(cal->GetNumberOfTimeBins());
292   AliGeomManager::LoadGeometry();
293
294
295   mgr->StartAnalysis("local",chain);
296
297   timer.Stop();
298   timer.Print();  
299
300   cal->Terminate();
301   delete field;
302   delete cdbManager;
303   for(Int_t it=NTRDTASKS-1; it>=0; it--) if(taskPtr[it]) delete taskPtr[it];
304   if(mcH) delete mcH;
305   delete esdH;
306   delete mgr;
307   delete chain;
308   if(MEM) delete mem;
309   if(MEM) TMemStatViewerGUI::ShowGUI();
310 }
311
312
313 TChain* CreateESDChain(const char* filename, Int_t nfiles)
314 {
315   // Create the chain
316   TChain* chain = new TChain("esdTree");
317
318   if(!filename){
319     chain->Add(Form("%s/AliESDs.root", gSystem->pwd()));
320     return chain;
321   }
322
323
324   // read ESD files from the input list.
325   ifstream in;
326   in.open(filename);
327   TString esdfile;
328   while(in.good() && (nfiles--) ) {
329     in >> esdfile;
330     if (!esdfile.Contains("root")) continue; // protection
331     chain->Add(esdfile.Data());
332   }
333
334   in.close();
335
336   return chain;
337 }