]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/macros/AddPerformanceTask.C
Update of ITS tracking check task and related macros
[u/mrichter/AliRoot.git] / PWG1 / macros / AddPerformanceTask.C
1 ///////////////////////////////////////////////////////////////////////////////
2 // Macro to setup AliPerformanceTask for 
3 // TPC performance to be run on QA train
4 //
5 // By default 6 performance components are added to 
6 // the output: 
7 // 1. AliPerformanceRes (TPC track resolution at DCA)
8 // 2. AliPerformanceResTPCInner (TPC track resolution at inner TPC wall)
9 // 3. AliPerformanceEff (TPC track reconstruction efficieny)
10 // 4. AliPerformanceDEdxTPCInner (TPC dEdx response - track parameters at TPC inner wall)
11 // 5. AliPerformanceDCA (TPC impact parameter resolution at DCA)
12 // 6. AliPerformanceTPC (TPC cluster and track information)
13 //
14 //
15 // To use these components one needs to have access to ESD and MC events 
16 // and track references.
17 //
18 // Usage on the analysis train:
19 // #include "AddPerformanceTask.h"
20 // #include "PWG1/macros/AddPerformanceTask.C"
21 //
22 // gSystem->Load("libPWG1.so");
23 //
24 // AliAnalysisManager *mgr = new AliAnalysisManager("Post Reconstruction QA");
25 // gROOT->LoadMacro("$ALICE_ROOT/PWG1/macros/AddPerformanceTask.C")
26 // AddPerformanceTask(mgr,"ALL");
27 // 
28 // 
29 // Output:
30 // TPC.PerformanceTask.root file with performance components is created.
31 //
32 // Each of the components contains THnSparse generic histograms which 
33 // have to be analysed by using Analyse() function. Each component contains such function.
34 //
35 //24.04.2009 -  J.Otwinowski@gsi.de
36 ///////////////////////////////////////////////////////////////////////////////
37
38 #if ! defined (__CINT__) || defined (__MAKECINT__)
39 #include <Riostream.h>
40
41 #include "TROOT.h"
42 #include "TClass.h"
43 #include "TError.h"
44
45 #include "AliLog.h"
46 #include "AliAnalysisManager.h"
47 #include "PWG1/AliPerformanceTask.h"
48 #include "PWG1/AliPerformanceObject.h"
49 #include "PWG1/AliPerformanceEff.h"
50 #include "PWG1/AliPerformanceDEdx.h"
51 #include "PWG1/AliPerformanceTPC.h"
52 #include "PWG1/AliPerformanceDCA.h"
53 #include "PWG1/AliPerformanceRes.h"
54 #include "PWG1/AliMCInfoCuts.h"
55 #include "PWG1/AliRecInfoCuts.h"
56 #include "PWG1/macros/AddPerformanceTask.h"
57 #endif
58
59 Int_t ParseTPC(Char_t *tpc);
60
61 //____________________________________________
62 void AddPerformanceTask(AliAnalysisManager *mgr=0, Char_t *tpc="ALL")
63 {
64   if(!mgr) { 
65     Error("AddPerformanceTask","AliAnalysisManager not set!");
66     return;
67   }
68   // parse options
69   Int_t fSteerTPC = ParseTPC(tpc);
70
71   // Add task
72   AliPerformanceTask *task = new AliPerformanceTask("Performance","TPC Performance");
73   if (!task) return;
74   if ( mgr->GetMCtruthEventHandler() ) task->SetUseMCInfo(kTRUE);
75   mgr->AddTask(task);
76
77   // Create TPC-ESD track reconstruction cuts
78   AliRecInfoCuts *pRecInfoCuts = new AliRecInfoCuts(); 
79   if(pRecInfoCuts) {
80     pRecInfoCuts->SetMaxDCAToVertexXY(3.0);
81     pRecInfoCuts->SetMaxDCAToVertexZ(3.0);
82     pRecInfoCuts->SetMinNClustersTPC(50);
83     pRecInfoCuts->SetMinNClustersITS(2);
84     pRecInfoCuts->SetHistogramsOn(kFALSE); 
85   }
86   // Create TPC-MC track reconstruction cuts
87   AliMCInfoCuts  *pMCInfoCuts = new AliMCInfoCuts();
88   if(pMCInfoCuts) {
89     pMCInfoCuts->SetMinTrackLength(50);
90   }
91
92   //
93   // Create performance objects for TPC and set cuts 
94   //
95
96   // TPC at DCA
97   TClass ctask; AliPerformanceObject *perf = 0x0;
98   for(Int_t icomp=0; icomp<NTPCTASKS; icomp++) {
99     if(!(TSTTPCBIT(fSteerTPC, icomp))) continue;
100     TString  s(fgkTPCtaskClassName[icomp]);
101     if(s.CompareTo("AliPerformanceEff") == 0) {
102       task->AddPerformanceObject((perf = new AliPerformanceEff(fgkTPCtaskClassName[icomp],fgkTPCtaskClassName[icomp],kTPCMode,fHpt)));
103     } else if (s.CompareTo("AliPerformanceRes") == 0) {
104       task->AddPerformanceObject((perf = new AliPerformanceRes(fgkTPCtaskClassName[icomp],fgkTPCtaskClassName[icomp],kTPCMode,fHpt)));
105     } else if (s.CompareTo("AliPerformanceTPC") == 0) {
106       task->AddPerformanceObject((perf = new AliPerformanceTPC(fgkTPCtaskClassName[icomp],fgkTPCtaskClassName[icomp],kTPCMode,fHpt)));
107     } else if (s.CompareTo("AliPerformanceDCA") == 0) {
108       task->AddPerformanceObject((perf = new AliPerformanceDCA(fgkTPCtaskClassName[icomp],fgkTPCtaskClassName[icomp],kTPCMode,fHpt)));
109     } else {
110       Warning("No TPC at DCA mode for ",fgkTPCtaskClassName[icomp]);
111     }
112     perf->SetAliMCInfoCuts(pMCInfoCuts);
113     perf->SetAliRecInfoCuts(pRecInfoCuts);
114   }
115
116   // TPC at inner wall
117   for(Int_t icomp=0; icomp<NTPCTASKS; icomp++) {
118     if(!(TSTTPCBIT(fSteerTPC, icomp))) continue;
119     TString  s(fgkTPCtaskClassName[icomp]);
120     if(s.CompareTo("AliPerformanceRes") == 0) {
121       task->AddPerformanceObject((perf = new AliPerformanceRes("AliPerformanceResTPCInner","AliPerformanceResTPCInner",kTPCInnerMode,fHpt)));
122     } else if (s.CompareTo("AliPerformanceDEdx") == 0) {
123       task->AddPerformanceObject((perf = new AliPerformanceDEdx("AliPerformanceDEdxTPCInner","AliPerformanceDEdxTPCInner",kTPCInnerMode,fHpt)));
124     } else {
125       Warning("No TPC at inner wall mode for ",fgkTPCtaskClassName[icomp]);
126     } 
127     perf->SetAliMCInfoCuts(pMCInfoCuts);
128     perf->SetAliRecInfoCuts(pRecInfoCuts);
129   }
130
131   // Create containers for input
132   AliAnalysisDataContainer *cinput = mgr->GetCommonInputContainer();
133   mgr->ConnectInput(task, 0, cinput);
134
135   // Create containers for output
136   AliAnalysisDataContainer *coutput = mgr->CreateContainer("coutput", TList::Class(), AliAnalysisManager::kOutputContainer, Form("TPC.%s.root", task->GetName()));
137   mgr->ConnectOutput(task, 0, coutput);
138
139   // Enable debug printouts
140   mgr->SetDebugLevel(0);
141 }
142
143 //____________________________________________
144 Int_t ParseTPC(Char_t *tpc)
145 {
146   Printf("---------------------------------------");
147   Printf("TPC Performance Task Configuration Options:");
148   Printf("---------------------------------------");
149   for(int i=0; i<NTPCTASKS+1; i++) {
150     Printf("%s",fgkTPCtaskOpt[i]);
151   } 
152   Printf("%s","HPT");
153   Printf("---------------------------------------");
154   Printf("Used Options:");
155   Printf("---------------------------------------");
156
157   Int_t fSteerTask = 0;
158   TObjArray *tasksArray = TString(tpc).Tokenize(" ");
159   for(Int_t isel = 0; isel < tasksArray->GetEntriesFast(); isel++){
160     TString s = (dynamic_cast<TObjString *>(tasksArray->UncheckedAt(isel)))->String();
161     if(s.CompareTo("ALL") == 0) {
162       Printf("%s","ALL");
163       for(Int_t itask = 0; itask < NTPCPERFORMANCE; itask++) SETTPCBIT(fSteerTask, itask);
164       continue;
165     }
166     else if(s.CompareTo("HPT") == 0) {
167       fHpt = kTRUE;
168       Printf("%s","HPT");
169     } else { 
170       Bool_t foundOpt = kFALSE;  
171       for(Int_t itask = 0; itask < NTPCTASKS; itask++) {
172         if(s.CompareTo(fgkTPCtaskOpt[itask]) != 0) continue;
173         SETTPCBIT(fSteerTask, itask); //SETTPCBIT(fSteerTask, 0);
174         foundOpt = kTRUE;
175         Printf("%s", fgkTPCtaskOpt[itask]);
176         break;
177       }
178       if(!foundOpt) Info("run.C", Form("TPC task %s not implemented (yet).", s.Data()));
179     }
180   }
181   Printf("---------------------------------------");
182   return fSteerTask;
183 }