]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/run.C
Add resolution class
[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 //     "PID"  : TRD PID - pion efficiency 
10 //     "PIDR" : TRD PID - reference data
11 // 
12 // Authors:
13 //   Alex Bercuci (A.Bercuci@gsi.de) 
14 //   Markus Fasel (m.Fasel@gsi.de) 
15
16 #define BIT(n)       (1 << (n))
17 #define SETBIT(n,i)  ((n) |= BIT(i))
18 #define TESTBIT(n,i) ((Bool_t)(((n) & BIT(i)) != 0))
19
20 const Int_t fknTasks = 3;
21 Char_t *fTaskName[fknTasks] = {"Barrel Tracking Effiency", "Combined Tracking Efficiency", "Tracking Resolution"};
22 enum AliTRDrecoTasks{
23   kTrackingEfficiency = 0
24   ,kTrackingCombinedEfficiency = 1
25   ,kTrackingResolution = 2
26 };
27 void run(Char_t *tasks="ALL", const Char_t *files=0x0, Int_t nmax=-1)
28 {
29
30
31   TStopwatch timer;
32   timer.Start();
33
34   gSystem->Load("libANALYSIS.so");
35   gSystem->Load("libTRDqaRec.so");
36   
37   Int_t fSteerTask = 0; 
38   TObjArray *task = TString(tasks).Tokenize(" ");
39   for(Int_t isel = 0; isel < task->GetEntriesFast(); isel++){
40     TString s = (dynamic_cast<TObjString *>(task->UncheckedAt(isel)))->String();
41     if(s.CompareTo("ALL") == 0){
42       for(Int_t itask = 0; itask < fknTasks; itask++) SETBIT(fSteerTask, itask);
43       continue;
44     } else if(s.CompareTo("EFF") == 0){
45       SETBIT(fSteerTask, kTrackingEfficiency);
46       continue;
47     } else if(s.CompareTo("EFFC") == 0){
48       SETBIT(fSteerTask, kTrackingCombinedEfficiency);
49       continue;
50     } else if(s.CompareTo("RES" ) == 0){
51       SETBIT(fSteerTask, kTrackingResolution);
52       continue;
53     } else{
54       Info("run.C", Form("Task %s not implemented (yet).", s.Data()));
55       continue;
56     }
57   }
58   printf("\n\tRUNNING TRAIN FOR TASKS:\n");
59   for(itask = 0; itask < fknTasks; itask++){
60     if(TESTBIT(fSteerTask, itask)) printf("\t%s\n", fTaskName[itask]);
61   }
62
63   //____________________________________________//
64   gROOT->LoadMacro(Form("%s/TRD/qaRec/CreateESDChain.C", gSystem->ExpandPathName("$ALICE_ROOT")));
65   TChain *chain = CreateESDChain(files, nmax);
66   //chain->SetBranchStatus("*", 0);
67   chain->SetBranchStatus("*FMD*",0);
68   chain->SetBranchStatus("*Calo*",0);
69   chain->SetBranchStatus("Tracks", 1);
70   chain->SetBranchStatus("ESDfriend*",1);
71   chain->Lookup();
72   chain->GetListOfFiles()->Print();
73   printf("\n ----> CHAIN HAS %d ENTRIES <----\n\n", (Int_t)chain->GetEntries());
74   
75   AliLog::SetGlobalLogLevel(AliLog::kError);
76
77   //____________________________________________
78   // Make the analysis manager
79   AliAnalysisManager *mgr = new AliAnalysisManager("TRD Track Info Manager");
80   //mgr->SetSpecialOutputLocation(source); // To Be Changed
81   AliVEventHandler* esdH = new AliESDInputHandler;
82   mgr->SetInputEventHandler(esdH);  
83   AliMCEventHandler *mc = new AliMCEventHandler();
84   mgr->SetMCtruthEventHandler(mc);
85   //mgr->SetDebugLevel(10);
86
87   //____________________________________________
88   // TRD track summary generator
89   AliTRDtrackInfoGen *task1 = new AliTRDtrackInfoGen();
90   task1->SetDebugLevel(1);
91   mgr->AddTask(task1);
92   // Create containers for input/output
93   AliAnalysisDataContainer *cinput1 = mgr->CreateContainer("data", TChain::Class(), AliAnalysisManager::kInputContainer);
94   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("TrackInfoList", TObjArray::Class(), AliAnalysisManager::kExchangeContainer);
95   mgr->ConnectInput( task1, 0, cinput1);
96   mgr->ConnectOutput(task1, 0, coutput1);
97
98   //____________________________________________
99   // TRD barrel tracking efficiency
100   if(TESTBIT(fSteerTask, kTrackingEfficiency)){
101     AliTRDtrackingEfficiency *task2 = new AliTRDtrackingEfficiency();
102     task2->SetDebugLevel(1);
103     mgr->AddTask(task2);
104     //Create containers for input/output
105     AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("TrackingEfficiency", TList::Class(), AliAnalysisManager::kOutputContainer, "TRD.TrackingEfficiency.root");
106     mgr->ConnectInput( task2, 0, coutput1);
107     mgr->ConnectOutput(task2, 0, coutput2);
108   }
109
110   //____________________________________________
111   // TRD combined tracking efficiency
112   if(TESTBIT(fSteerTask, kTrackingCombinedEfficiency)){
113     AliTRDtrackingEfficiencyCombined *task3 = new AliTRDtrackingEfficiencyCombined();
114     task3->SetDebugLevel(0);
115     mgr->AddTask(task3);
116     // Create containers for input/output
117     AliAnalysisDataContainer *coutput3 = mgr->CreateContainer("TrackingEfficiencyCombined", TObjArray::Class(), AliAnalysisManager::kOutputContainer, "TRD.TrackingEfficiencyCombined.root");
118     mgr->ConnectInput( task3, 0, coutput1);
119     mgr->ConnectOutput(task3, 0, coutput3);
120   }
121
122   //____________________________________________
123   // TRD combined tracking efficiency
124   if(TESTBIT(fSteerTask, kTrackingResolution)){
125     AliTRDtrackingResolution *task4 = new AliTRDtrackingResolution();
126     task4->SetDebugLevel(1);
127     mgr->AddTask(task4);
128     // Create containers for input/output
129     AliAnalysisDataContainer *coutput4 = mgr->CreateContainer("Tracking Resolution", TList::Class(), AliAnalysisManager::kOutputContainer, "TRD.TrackingResolution.root");
130     mgr->ConnectInput( task4, 0, coutput1);
131     mgr->ConnectOutput(task4, 0, coutput4);
132   }
133
134   if (!mgr->InitAnalysis()) return;
135   mgr->PrintStatus();
136   mgr->StartAnalysis("local",chain);
137
138   timer.Stop();
139   timer.Print();
140 }