]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/run.C
update train
[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)
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 //     "ALGN" : TRD alignment
12 //     "PID"  : TRD PID - pion efficiency 
13 //     "PIDR" : TRD PID - reference data
14 //     "DET"  : Basic TRD Detector checks
15 //     "NOFR" : Data set does not have AliESDfriends.root 
16 //     "NOMC" : Data set does not have Monte Carlo Informations (real data), so all tasks which rely
17 //              on MC information are switched off
18 //
19 // In compiled mode : 
20 // Don't forget to load first the libraries
21 // gSystem->Load("libMemStat.so")
22 // gSystem->Load("libMemStatGui.so")
23 // gSystem->Load("libANALYSIS.so")
24 // gSystem->Load("libTRDqaRec.so")
25 // gSystem->Load("libNetx.so") ;
26 // gSystem->Load("libRAliEn.so");
27 //
28 // Authors:
29 //   Alex Bercuci (A.Bercuci@gsi.de) 
30 //   Markus Fasel (m.Fasel@gsi.de) 
31
32 #ifndef __CINT__
33 #include <Riostream.h>
34
35 #include "TStopwatch.h"
36 #include "TMemStat.h"
37 #include "TMemStatViewerGUI.h"
38
39 #include "TROOT.h"
40 #include "TSystem.h"
41 #include "TError.h"
42 #include "TChain.h"
43 #include "TGrid.h"
44 #include "TAlienCollection.h"
45 #include "TGridCollection.h"
46 #include "TGridResult.h"
47
48 #include "AliMagFMaps.h"
49 #include "AliTracker.h"
50 #include "AliLog.h"
51 #include "AliCDBManager.h"
52 #include "AliGeomManager.h"
53 #include "AliAnalysisManager.h"
54 #include "AliAnalysisDataContainer.h"
55 #include "AliMCEventHandler.h"
56 #include "AliESDInputHandler.h"
57
58 #include "TRD/AliTRDtrackerV1.h"
59 #include "TRD/AliTRDcalibDB.h"
60 #include "TRD/qaRec/AliTRDtrackInfo/AliTRDeventInfo.h"
61 #include "TRD/qaRec/AliTRDtrackInfoGen.h"
62 #include "TRD/qaRec/AliTRDtrackingEfficiency.h"
63 #include "TRD/qaRec/AliTRDtrackingEfficiencyCombined.h"
64 #include "TRD/qaRec/AliTRDtrackingResolution.h"
65 #include "TRD/qaRec/AliTRDcalibration.h"
66 #include "TRD/qaRec/AliTRDalignmentTask.h"
67 #include "TRD/qaRec/AliTRDpidChecker.h"
68 #include "TRD/qaRec/AliTRDpidRefMaker.h"
69 #include "TRD/qaRec/AliTRDcheckDetector.h"
70 #include "TRD/qaRec/AliTRDclusterResolution.h"
71 #endif
72
73 #include "run.h"
74
75 Bool_t MEM = kFALSE;
76
77 TChain* MakeChainLST(const char* filename = 0x0);
78 TChain* MakeChainXML(const char* filename = 0x0);
79 void run(Char_t *tasks="ALL", const Char_t *files=0x0)
80 {
81   TMemStat *mem = 0x0;
82   if(MEM){ 
83     gSystem->Load("libMemStat.so");
84     gSystem->Load("libMemStatGui.so");
85     mem = new TMemStat("new, gnubuildin");
86     mem->AddStamp("Start");
87   }
88
89   TStopwatch timer;
90   timer.Start();
91
92   if(gSystem->Load("libANALYSIS.so")<0) return;
93   if(gSystem->Load("libTRDqaRec.so")<0) return;
94   
95   Bool_t fHasMCdata = kTRUE;
96   Bool_t fHasFriends = kTRUE;
97   TObjArray *tasksArray = TString(tasks).Tokenize(" ");
98
99   Int_t fSteerTask = 0; SETBIT(fSteerTask, kInfoGen);
100   for(Int_t isel = 0; isel < tasksArray->GetEntriesFast(); isel++){
101     TString s = (dynamic_cast<TObjString *>(tasksArray->UncheckedAt(isel)))->String();
102     if(s.CompareTo("ALL") == 0){
103       for(Int_t itask = 1; itask < NQATASKS; itask++) SETBIT(fSteerTask, itask);
104       continue;
105     } else if(s.CompareTo("NOFR") == 0){ 
106       fHasFriends = kFALSE;
107     } else if(s.CompareTo("NOMC") == 0){ 
108       fHasMCdata = kFALSE;
109     } else { 
110       Bool_t foundOpt = kFALSE;  
111       for(Int_t itask = 1; itask < NTRDTASKS; itask++){
112         if(s.CompareTo(fgkTRDtaskOpt[itask]) != 0) continue;
113         SETBIT(fSteerTask, itask);
114         foundOpt = kTRUE;
115         break;
116       }
117       if(!foundOpt) Info("run.C", Form("Task %s not implemented (yet).", s.Data()));
118     }
119   }
120   // extra rules for calibration tasks
121   if(TSTBIT(fSteerTask, kClErrParam)) SETBIT(fSteerTask, kResolution);
122   if(TSTBIT(fSteerTask, kPIDRefMaker)) SETBIT(fSteerTask, kPIDChecker);
123   if(TSTBIT(fSteerTask, kAlignment)) SETBIT(fSteerTask, kResolution);
124
125   // define task list pointers;
126   AliTRDrecoTask *taskPtr[NTRDTASKS], *task = 0x0;
127   memset(taskPtr, 0, NTRDTASKS*sizeof(AliAnalysisTask*));
128
129   //____________________________________________//
130   // DEFINE DATA CHAIN
131   TChain *chain = 0x0;
132   if(!files) chain = MakeChainLST();
133   else{
134     TString fn(files);
135     if(fn.EndsWith("xml")) chain = MakeChainXML(files);
136     else chain = MakeChainLST(files);
137   }
138   if(!chain) return;
139
140   //chain->SetBranchStatus("*", 0);
141   chain->SetBranchStatus("*FMD*",0);
142   chain->SetBranchStatus("*Calo*",0);
143   chain->SetBranchStatus("Tracks", 1);
144   chain->SetBranchStatus("ESDfriend*",1);
145   chain->Lookup();
146   chain->GetListOfFiles()->Print();
147   printf("\n ----> CHAIN HAS %d ENTRIES <----\n\n", (Int_t)chain->GetEntries());
148   
149   AliLog::SetGlobalLogLevel(AliLog::kError);
150
151   //____________________________________________
152   // Make the analysis manager
153   AliAnalysisManager *mgr = new AliAnalysisManager("TRD QA Reconstruction Manager");
154   //mgr->SetSpecialOutputLocation(source); // To Be Changed
155   AliVEventHandler *esdH = 0x0, *mcH = 0x0;
156   mgr->SetInputEventHandler(esdH = new AliESDInputHandler);
157   if(fHasMCdata) mgr->SetMCtruthEventHandler(mcH = new AliMCEventHandler());
158   //mgr->SetDebugLevel(10);
159
160   //____________________________________________
161   // TRD track summary generator
162   mgr->AddTask(task = new AliTRDtrackInfoGen());
163   taskPtr[(Int_t)kInfoGen] = task;
164   task->SetDebugLevel(1);
165   task->SetMCdata(fHasMCdata);
166   // Create containers for input/output
167   AliAnalysisDataContainer *cinput1 = mgr->CreateContainer("data", TChain::Class(), AliAnalysisManager::kInputContainer);
168   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("trackInfo", TObjArray::Class(), AliAnalysisManager::kExchangeContainer);
169   AliAnalysisDataContainer *coutput1a = mgr->CreateContainer("eventInfo", AliTRDeventInfo::Class(), AliAnalysisManager::kExchangeContainer);
170   mgr->ConnectInput( task, 0, cinput1);
171   mgr->ConnectOutput(task, 0, coutput1);
172   mgr->ConnectOutput(task, 1, coutput1a);
173
174   //____________________________________________
175   // TRD detector checker
176         if(TSTBIT(fSteerTask, kCheckDetector)){
177     mgr->AddTask(task = new AliTRDcheckDetector());
178     taskPtr[(Int_t)kCheckDetector] = task;
179     task->SetDebugLevel(4);
180     task->SetMCdata(fHasMCdata);
181     
182     // Create containers for input/output
183     mgr->ConnectInput( task, 0, coutput1);
184     mgr->ConnectInput( task, 1, coutput1a);
185     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
186   }
187
188   //____________________________________________
189   // TRD barrel tracking efficiency
190   if(fHasMCdata && TSTBIT(fSteerTask, kTrackingEff)){
191     mgr->AddTask(task = new AliTRDtrackingEfficiency());
192     taskPtr[(Int_t)kTrackingEff] = task;
193     task->SetDebugLevel(0);
194
195     //Create containers for input/output
196     mgr->ConnectInput( task, 0, coutput1);
197     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
198   }
199
200   //____________________________________________
201   // TRD combined tracking efficiency
202   if(fHasMCdata && TSTBIT(fSteerTask, kTrackingEffMC)){
203     mgr->AddTask(task = new AliTRDtrackingEfficiencyCombined());
204     taskPtr[(Int_t)kTrackingEffMC] = task;
205     task->SetDebugLevel(0);
206
207     // Create containers for input/output
208     mgr->ConnectInput( task, 0, coutput1);
209     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
210   }
211
212   //____________________________________________
213   // TRD tracking resolution
214   if(TSTBIT(fSteerTask, kResolution)){
215     mgr->AddTask(task = new AliTRDtrackingResolution());
216     taskPtr[(Int_t)kResolution] = task;
217     task->SetMCdata(fHasMCdata);
218     task->SetPostProcess(kFALSE);
219     task->SetDebugLevel(1);
220     
221     // Create containers for input/output
222     mgr->ConnectInput( task, 0, coutput1);
223     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
224
225     // Create output containers for calibration tasks
226     const Int_t nc = 5;
227     const Char_t *cn[nc] = {"Cl", "TrkltY", "TrkltPhi", "ClMC", "TrkltYMC"}; 
228     AliAnalysisDataContainer *co[nc]; 
229     for(Int_t ic = 0; ic<nc; ic++){
230       co[ic] = mgr->CreateContainer(Form("%s%s", task->GetName(), cn[ic]), TObjArray::Class(), AliAnalysisManager::kExchangeContainer);
231       mgr->ConnectOutput(task, 1+ic, co[ic]);
232     }
233     
234     // test reconstruction calibration plugin
235     if(TSTBIT(fSteerTask, kClErrParam)){
236       mgr->AddTask(task = new AliTRDclusterResolution());
237       taskPtr[(Int_t)kClErrParam] = task;
238       mgr->ConnectInput(task, 0, co[0]);
239       mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
240   
241       mgr->AddTask(task = new AliTRDclusterResolution());
242       taskPtr[(Int_t)kClErrParam+1] = task;
243       mgr->ConnectInput(task, 0, co[3]);
244       mgr->ConnectOutput(task, 0, mgr->CreateContainer(Form("%sMC", task->GetName()), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sMC.root", task->GetName())));
245     }
246   }
247
248   //____________________________________________
249   // TRD calibration
250   if(TSTBIT(fSteerTask, kCalibration)){
251     mgr->AddTask(task = new AliTRDcalibration());
252     taskPtr[(Int_t)kCalibration] = task;
253     ((AliTRDcalibration*)task)->SetLow(0);
254     ((AliTRDcalibration*)task)->SetHigh(30);
255     ((AliTRDcalibration*)task)->SetFillZero(kFALSE);
256     task->SetDebugLevel(0);
257
258     // Create containers for input/output
259     mgr->ConnectInput(task, 0, coutput1);
260     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
261   }
262   
263   //____________________________________________
264   // TRD alignment
265   if(TSTBIT(fSteerTask, kAlignment)){
266     mgr->AddTask(task = new AliTRDalignmentTask());
267     taskPtr[(Int_t)kAlignment] = task;
268     task->SetDebugLevel(0);
269
270     // Create containers for input/output
271     mgr->ConnectInput(task, 0, coutput1);
272     mgr->ConnectOutput(task, 0, mgr->CreateContainer(Form("h%s", task->GetName()), TObjArray::Class(), AliAnalysisManager::kExchangeContainer));
273
274     mgr->ConnectOutput(task, 1, mgr->CreateContainer(task->GetName(), TTree::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
275   }
276   
277   //____________________________________________
278   // TRD PID
279   if(TSTBIT(fSteerTask, kPIDChecker)){
280     mgr->AddTask(task = new AliTRDpidChecker());
281     taskPtr[(Int_t)kPIDChecker] = task;
282     task->SetDebugLevel(0);
283     task->SetMCdata(fHasMCdata);
284     
285     // Create containers for input/output
286     mgr->ConnectInput( task, 0, coutput1);
287     mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
288
289     //____________________________________________
290     // TRD pid reference 
291     if(TSTBIT(fSteerTask, kPIDRefMaker)){
292       mgr->AddTask(task = new AliTRDpidRefMaker());
293       taskPtr[(Int_t)kPIDRefMaker] = task;
294       task->SetDebugLevel(0);
295       task->SetMCdata(fHasMCdata);
296       
297       // Create containers for input/output
298       mgr->ConnectInput( task, 0, coutput1);
299       mgr->ConnectOutput(task, 0, mgr->CreateContainer(task->GetName(), TObjArray::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%s.root", task->GetName())));
300       mgr->ConnectOutput(task, 1, mgr->CreateContainer(Form("%sNN", task->GetName()), TTree::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sNN.root", task->GetName())));
301       mgr->ConnectOutput(task, 2, mgr->CreateContainer(Form("%sLQ", task->GetName()), TTree::Class(), AliAnalysisManager::kOutputContainer, Form("TRD.Task%sLQ.root", task->GetName())));
302     }
303   }
304
305
306   if (!mgr->InitAnalysis()) return;
307   printf("\n\tRUNNING TRAIN FOR TASKS:\n");
308   for(Int_t itask = 1; itask < NTRDTASKS; itask++){
309     if(TSTBIT(fSteerTask, itask)) printf("\t   %s [%s]\n",  taskPtr[itask]->GetName(), taskPtr[itask]->GetTitle());
310   }
311   printf("\n\n");
312   //mgr->PrintStatus();
313
314
315   AliCDBManager *cdbManager = AliCDBManager::Instance();
316   cdbManager->SetDefaultStorage("local://$ALICE_ROOT");
317   cdbManager->SetRun(0);
318   cdbManager->SetCacheFlag(kFALSE);
319  
320   // initialize TRD settings
321   AliMagFMaps *field = new AliMagFMaps("Maps","Maps", 2, 1., 10., AliMagFMaps::k5kG);
322   AliTracker::SetFieldMap(field, kTRUE);
323   AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
324   AliTRDtrackerV1::SetNTimeBins(cal->GetNumberOfTimeBins());
325   AliGeomManager::LoadGeometry();
326
327
328   mgr->StartAnalysis("local",chain);
329
330   timer.Stop();
331   timer.Print();  
332
333   cal->Terminate();
334   delete field;
335   delete cdbManager;
336   for(Int_t it=NTRDTASKS; it--; ){ 
337     if(taskPtr[it]){ 
338       printf("Cleaning %s [%s] ...\n", fgkTRDtaskClassName[it], taskPtr[it]->GetTitle());
339       delete taskPtr[it];
340     }
341   }
342   if(mcH) delete mcH;
343   delete esdH;
344   delete mgr;
345   delete chain;
346   if(MEM) delete mem;
347   if(MEM) TMemStatViewerGUI::ShowGUI();
348 }
349
350 //____________________________________________
351 TChain* MakeChainLST(const char* filename)
352 {
353   // Create the chain
354   TChain* chain = new TChain("esdTree");
355
356   if(!filename){
357     chain->Add(Form("%s/AliESDs.root", gSystem->pwd()));
358     return chain;
359   }
360
361
362   // read ESD files from the input list.
363   ifstream in;
364   in.open(filename);
365   TString esdfile;
366   while(in.good()) {
367     in >> esdfile;
368     if (!esdfile.Contains("root")) continue; // protection
369     chain->Add(esdfile.Data());
370   }
371
372   in.close();
373
374   return chain;
375 }
376
377 //____________________________________________
378 TChain* MakeChainXML(const char* xmlfile)
379 {
380   if (!TFile::Open(xmlfile)) {
381     Error("MakeChainXML", Form("No file %s was found", xmlfile));
382     return 0x0;
383   }
384
385   if(gSystem->Load("libNetx.so")<0) return 0x0;
386   if(gSystem->Load("libRAliEn.so")<0) return 0x0;
387   TGrid::Connect("alien://") ;
388
389   TGridCollection *collection = (TGridCollection*) TAlienCollection::Open(xmlfile);
390   if (!collection) {
391     Error("MakeChainXML", Form("No collection found in %s", xmlfile)) ; 
392     return 0x0; 
393   }
394   //collection->CheckIfOnline();
395
396   TGridResult* result = collection->GetGridResult("",0 ,0);
397   if(!result->GetEntries()){
398     Error("MakeChainXML", Form("No entries found in %s", xmlfile)) ; 
399     return 0x0; 
400   }
401   // Makes the ESD chain 
402   TChain* chain = new TChain("esdTree");
403   for (Int_t idx = 0; idx < result->GetEntries(); idx++) {
404     chain->Add(result->GetKey(idx, "turl")); 
405   }
406   return chain;
407 }