]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliComparisonTask.cxx
Removed hidden symbols (Marian)
[u/mrichter/AliRoot.git] / PWG1 / AliComparisonTask.cxx
1 //------------------------------------------------------------------------------\r
2 // Implementation of the AliComparisonTask class. It compares properties of the \r
3 // reconstructed and MC particle tracks under several conditions. \r
4 // As the input it requires the TTree with AliRecInfo and AliMCInfo branches. \r
5 // The comparison output histograms are stored \r
6 // in the comparison objects: AliComparisonRes, AliComparisonEff, \r
7 // AliComparisonDEdx and AliComparisonDCA. Each of these objects also contains \r
8 // selection cuts which were used during filling the histograms.\r
9 // \r
10 // Author: J.Otwinowski 04/02/2008 \r
11 //------------------------------------------------------------------------------\r
12 \r
13 #include "iostream"\r
14 \r
15 #include "TChain.h"\r
16 #include "TTree.h"\r
17 #include "TH1F.h"\r
18 #include "TCanvas.h"\r
19 #include "TList.h"\r
20 #include "TFile.h"\r
21 \r
22 #include "AliAnalysisTask.h"\r
23 #include "AliAnalysisManager.h"\r
24 #include "AliESDEvent.h"\r
25 #include "AliESDInputHandler.h"\r
26 #include "AliESDVertex.h"\r
27 #include "AliMagFMaps.h"\r
28 #include "AliTracker.h"\r
29 #include "AliGeomManager.h"\r
30 \r
31 #include "AliMCInfo.h"\r
32 #include "AliESDRecInfo.h"\r
33 #include "AliMCInfoCuts.h"\r
34 #include "AliRecInfoCuts.h"\r
35 #include "AliComparisonRes.h"\r
36 #include "AliComparisonEff.h"\r
37 #include "AliComparisonDEdx.h"\r
38 #include "AliComparisonDCA.h"\r
39 #include "AliComparisonTask.h"\r
40 \r
41 using namespace std;\r
42 \r
43 ClassImp(AliComparisonTask)\r
44 \r
45 Int_t AliComparisonTask::evtNumber = 0;\r
46 \r
47 //_____________________________________________________________________________\r
48 AliComparisonTask::AliComparisonTask(const char *name) \r
49   : AliAnalysisTask(name, "")\r
50   , fTree(0)\r
51   , fInfoMC(0)\r
52   , fInfoRC(0)\r
53   , fCompRes(0)\r
54   , fCompEff(0)\r
55   , fCompDEdx(0)\r
56   , fCompDCA(0)\r
57   , fOutput(0)\r
58   , fMagField(0)\r
59   , fMagFMap(0)\r
60   , fGeom(0)\r
61 {\r
62   // Constructor\r
63 \r
64   // Define input and output slots here\r
65   DefineInput(0, TChain::Class());\r
66   DefineOutput(0, TList::Class());\r
67 \r
68   // set default mag. field\r
69   SetMagField();\r
70   \r
71   // set default geometry\r
72   SetGeometry();\r
73 }\r
74 \r
75 //_____________________________________________________________________________\r
76 AliComparisonTask::~AliComparisonTask()\r
77 {\r
78   if(fOutput)   delete fOutput;  fOutput =0; \r
79   if(fMagFMap)  delete fMagFMap;  fMagFMap =0; \r
80 }\r
81 \r
82 //_____________________________________________________________________________\r
83 void AliComparisonTask::ConnectInputData(Option_t *) \r
84 {\r
85   // Connect input data \r
86   // Called once\r
87 \r
88   fTree = dynamic_cast<TTree*> (GetInputData(0));\r
89   if (!fTree) {\r
90     Printf("ERROR: Could not read chain from input slot 0");\r
91   } else {\r
92     fTree->SetBranchStatus("*",1);\r
93   }\r
94 \r
95   if(fTree->GetBranch("MC") &&  fTree->GetBranch("RC")) {\r
96     fTree->GetBranch("MC")->SetAddress(&fInfoMC);\r
97     fTree->GetBranch("RC")->SetAddress(&fInfoRC);\r
98   } else {\r
99       Printf("ERROR: Could not get MC and RC branches");\r
100   }\r
101   \r
102   // set mag. field map \r
103   fMagFMap = new AliMagFMaps("Maps","Maps", 2, 1., 10., fMagField);\r
104   AliTracker::SetFieldMap(fMagFMap,kFALSE);\r
105 \r
106   // set geommetry\r
107   AliGeomManager::LoadGeometry(fGeom);\r
108 }\r
109 \r
110 //_____________________________________________________________________________\r
111 void AliComparisonTask::CreateOutputObjects()\r
112 {\r
113   // Create histograms\r
114   // Called once\r
115 \r
116   fOutput = new TList;\r
117   fOutput->SetOwner();\r
118 \r
119   if(fCompRes) fOutput->Add(fCompRes);\r
120   else \r
121      Printf("WARNING: AliComparisonRes is not added to the output");\r
122 \r
123   if(fCompEff) fOutput->Add(fCompEff);\r
124   else \r
125     Printf("WARNING: AliComparisonEff is not added to the output");\r
126 \r
127   if(fCompDEdx) fOutput->Add(fCompDEdx);\r
128   else \r
129     Printf("WARNING: AliComparisonDEdx is not added to the output");\r
130 \r
131   if(fCompDCA) fOutput->Add(fCompDCA);\r
132   else \r
133      Printf("WARNING: AliComparisonDCA is not added to the output");\r
134 }\r
135 \r
136 //_____________________________________________________________________________\r
137 Bool_t AliComparisonTask::ReadEntry(Int_t evt) \r
138 {\r
139   Long64_t centry = fTree->LoadTree(evt);\r
140   if(centry < 0) return kFALSE;\r
141 \r
142   if(fTree->GetBranch("MC") &&  fTree->GetBranch("RC")) {\r
143     fTree->GetBranch("MC")->SetAddress(&fInfoMC);\r
144     fTree->GetBranch("RC")->SetAddress(&fInfoRC);\r
145   } else {\r
146       Printf("ERROR: Could not get MC and RC branches");\r
147           return kFALSE;\r
148   }\r
149   fTree->GetEntry(evt);\r
150 \r
151 return kTRUE;\r
152 }\r
153 //_____________________________________________________________________________\r
154 void AliComparisonTask::Exec(Option_t *) \r
155 {\r
156   // Main loop\r
157   // Called for each event\r
158 \r
159   if (!fInfoMC && !fInfoRC) {\r
160     Printf("ERROR: fInfoMC && fInfoRC not available");\r
161     return;\r
162   }\r
163 \r
164   // Process comparison\r
165   Bool_t status = ReadEntry(evtNumber);\r
166   if(status == kTRUE) \r
167   {\r
168      if(fCompRes)  fCompRes->Exec(fInfoMC,fInfoRC);\r
169      if(fCompEff)  fCompEff->Exec(fInfoMC,fInfoRC);\r
170      if(fCompDEdx) fCompDEdx->Exec(fInfoMC,fInfoRC);\r
171      if(fCompDCA)  fCompDCA->Exec(fInfoMC,fInfoRC);\r
172   }\r
173 \r
174   if( !( evtNumber % 10000) ) { \r
175     cout << evtNumber << endl;\r
176   }\r
177   evtNumber++;\r
178 \r
179   // Post output data.\r
180   PostData(0, fOutput);\r
181 }\r
182 \r
183 //_____________________________________________________________________________\r
184 void AliComparisonTask::Terminate(Option_t *) \r
185 {\r
186   // Called once at the end of the event loop\r
187   cout << "Terminate " << endl;\r
188 \r
189   TFile *out = new TFile("Output.root","RECREATE");\r
190   out->cd();\r
191 \r
192   fOutput = dynamic_cast<TList*> (GetOutputData(0));\r
193   if (!fOutput) {\r
194     Printf("ERROR: fOutput not available");\r
195     return;\r
196   }\r
197 \r
198   fCompRes = dynamic_cast<AliComparisonRes*> (fOutput->FindObject("AliComparisonRes"));\r
199   if (!fCompRes) {\r
200     Printf("WARNING: AliComparisonRes not available");\r
201     return;\r
202   }\r
203 \r
204   fCompEff = dynamic_cast<AliComparisonEff*> (fOutput->FindObject("AliComparisonEff"));\r
205   if (!fCompEff) {\r
206     Printf("WARNING: AliComparisonEff not available");\r
207     return;\r
208   }\r
209    \r
210   fCompDEdx = dynamic_cast<AliComparisonDEdx*> (fOutput->FindObject("AliComparisonDEdx"));\r
211   if (!fCompDEdx) {\r
212     Printf("WARNING: AliComparisonDEdx not available");\r
213     return;\r
214   }\r
215 \r
216   fCompDCA = dynamic_cast<AliComparisonDCA*> (fOutput->FindObject("AliComparisonDCA"));\r
217   if (!fCompDCA) {\r
218     Printf("WARNING: AliComparisonDCA not available");\r
219     return;\r
220   }\r
221 \r
222   fOutput->Write();\r
223   out->Close();\r
224 }\r