]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDrecoTask.cxx
updates for running on grid (Markus)
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDrecoTask.cxx
1 ///////////////////////////////////////////////////////
2 //
3 // Basic class for Performance/Calibration TRD tasks
4 // 
5 // It performs generic tasks like :
6 //   - data file manegment
7 //   - reference container management
8 //   - debug container management
9 //   - interaction with AliAnalysisManager
10 //   - Plot functor loop
11 //
12 // Author: Alexandru Bercuci <A.Bercuci@gsi.de>, 10/09/2008
13 //
14 /////////////////////////////////////////////////////////
15
16 #include "TClass.h"
17 #include "TMethod.h"
18 #include "TMethodCall.h"
19 #include "TMethodArg.h"
20 #include "TFile.h"
21 #include "TChain.h"
22 #include "TList.h"
23 #include "TMap.h"
24 #include "TH1.h"
25 #include "TF1.h"
26 #include "TObjArray.h"
27 #include "TDirectory.h"
28 #include "TTreeStream.h"
29
30 #include "AliLog.h"
31 #include "AliAnalysisTask.h"
32
33 #include "AliTRDrecoTask.h"
34
35 ClassImp(AliTRDrecoTask)
36
37 TList* AliTRDrecoTask::fgTrendPoint(NULL);
38 TTreeSRedirector* AliTRDrecoTask::fgDebugStream(NULL);
39 //_______________________________________________________
40 AliTRDrecoTask::AliTRDrecoTask()
41   : AliAnalysisTaskSE()
42   ,fNRefFigures(0)
43   ,fContainer(NULL)
44   ,fTracks(NULL)
45   ,fkTrack(NULL)
46   ,fkMC(NULL)
47   ,fkESD(NULL)
48   ,fPlotFuncList(NULL)
49 {
50 // Default constructor  
51 }
52
53 //_______________________________________________________
54 AliTRDrecoTask::AliTRDrecoTask(const char *name, const char *title)
55   : AliAnalysisTaskSE(name)
56   ,fNRefFigures(0)
57   ,fContainer(NULL)
58   ,fTracks(NULL)
59   ,fkTrack(NULL)
60   ,fkMC(NULL)
61   ,fkESD(NULL)
62   ,fPlotFuncList(NULL)
63 {
64 // Constructor for all derived performance tasks
65
66   SetTitle(title);
67   DefineInput (1, TObjArray::Class()); // track list
68   DefineOutput(1, TObjArray::Class()); // histogram list
69 }
70
71 //_______________________________________________________
72 AliTRDrecoTask::~AliTRDrecoTask() 
73 {
74
75   // Generic task destructor
76
77   AliDebug(2, Form(" Ending task %s[%s]", GetName(), GetTitle()));
78   if(fgDebugStream){ 
79     delete fgDebugStream;
80     fgDebugStream = NULL;
81   }
82
83   if(fPlotFuncList){
84     fPlotFuncList->Delete();
85     delete fPlotFuncList;
86     fPlotFuncList = NULL;
87   }
88   
89   if(fContainer){
90     if(fContainer->IsOwner()) fContainer->Delete();
91     delete fContainer;
92     fContainer = NULL;
93   }
94
95   if(fgTrendPoint){
96     TFile::Open("TRD.PerformanceTrend.root", "UPDATE");
97     fgTrendPoint->Write();
98     delete fgTrendPoint;
99     fgTrendPoint=NULL;
100     gFile->Close();
101   }
102 }
103
104 //_______________________________________________________
105 void AliTRDrecoTask::ConnectInputData(Option_t *)
106 {
107   //
108   // Connect input data
109   //
110   AliAnalysisTaskSE::ConnectInputData();
111   fTracks = dynamic_cast<TObjArray *>(GetInputData(1));
112 }
113
114 //_______________________________________________________
115 void AliTRDrecoTask::UserExec(Option_t *)
116 {
117 // Loop over Plot functors published by particular tasks
118
119   GetInputData(1);
120   if(!fPlotFuncList){
121     AliWarning("No functor list defined for the reference plots");
122     return;
123   }
124   if(!fTracks) return;
125   if(!fTracks->GetEntriesFast()) return;
126   else AliDebug(2, Form("Tracks[%d] for %s", fTracks->GetEntriesFast(), GetName()));
127
128   AliTRDtrackInfo *trackInfo = NULL;
129   TIter plotIter(fPlotFuncList);
130   TObjArrayIter trackIter(fTracks);
131   while((trackInfo = dynamic_cast<AliTRDtrackInfo*>(trackIter()))){
132     fkTrack = trackInfo->GetTrack();
133     fkMC    = trackInfo->GetMCinfo();
134     fkESD   = trackInfo->GetESDinfo();
135
136     TMethodCall *plot = NULL;
137     plotIter.Reset();
138     while((plot=dynamic_cast<TMethodCall*>(plotIter()))){
139       plot->Execute(this);
140     }
141   }
142   PostData(1, fContainer);
143 }
144
145 //_______________________________________________________
146 Bool_t AliTRDrecoTask::GetRefFigure(Int_t /*ifig*/)
147 {
148   AliWarning("Retrieving reference figures not implemented.");
149   return kFALSE;
150 }
151
152 //_______________________________________________________
153 Bool_t AliTRDrecoTask::PutTrendValue(const Char_t *name, Double_t val)
154 {
155 // Generic publisher for trend values
156
157   if(!fgTrendPoint){
158     fgTrendPoint = new TList();
159     fgTrendPoint->SetOwner();
160   }
161   fgTrendPoint->AddLast(new TNamed(Form("%s_%s", GetName(), name), Form("%f", val)));
162   return kTRUE;
163 }
164
165 //_______________________________________________________
166 void AliTRDrecoTask::InitFunctorList()
167 {
168 // Initialize list of functors
169
170   TClass *c = this->IsA();
171
172   TMethod *m = NULL;
173   TIter methIter(c->GetListOfMethods());
174   while((m=dynamic_cast<TMethod*>(methIter()))){
175     TString name(m->GetName());
176     if(!name.BeginsWith("Plot")) continue;
177     if(!fPlotFuncList) fPlotFuncList = new TList();
178     fPlotFuncList->AddLast(new TMethodCall(c, (const char*)name, ""));
179   }
180 }
181
182 //_______________________________________________________
183 Bool_t AliTRDrecoTask::Load(const Char_t *filename)
184 {
185 // Generic container loader
186
187   if(!TFile::Open(filename)){
188     AliWarning(Form("Couldn't open file %s.", filename));
189     return kFALSE;
190   }
191   TObjArray *o = NULL;
192   if(!(o = (TObjArray*)gFile->Get(GetName()))){
193     AliWarning("Missing histogram container.");
194     return kFALSE;
195   }
196   fContainer = (TObjArray*)o->Clone(GetName());
197   gFile->Close();
198   return kTRUE;
199 }
200
201 //________________________________________________________
202 Bool_t AliTRDrecoTask::Save(TObjArray * const results){
203   //
204   // Store the output graphs in a ROOT file
205   // Input TObject array will not be written as Key to the file,
206   // only content itself
207   //
208
209   TDirectory *cwd = gDirectory;
210   if(!TFile::Open(Form("TRD.Result%s.root", GetName()), "RECREATE")) return kFALSE;
211
212   TIterator *iter = results->MakeIterator();
213   TObject *inObject = NULL, *outObject = NULL;
214   while((inObject = iter->Next())){
215     outObject = inObject->Clone();
216     outObject->Write(NULL, TObject::kSingleKey);
217   }
218   delete iter;
219   gFile->Close(); delete gFile;
220   cwd->cd(); 
221   return kTRUE;
222 }
223
224 //_______________________________________________________
225 Bool_t AliTRDrecoTask::PostProcess()
226 {
227 // To be implemented by particular tasks
228
229   AliWarning("Post processing of reference histograms not implemented.");
230   return kFALSE;
231 }
232
233 //_______________________________________________________
234 void AliTRDrecoTask::SetDebugLevel(Int_t level)
235 {
236 // Generic debug handler
237
238   AliAnalysisTaskSE::SetDebugLevel(level);
239   if(DebugLevel()>=1){
240     TDirectory *savedir = gDirectory;
241     fgDebugStream = new TTreeSRedirector("TRD.DebugPerformance.root");
242     savedir->cd();
243   }
244 }
245
246 //____________________________________________________________________
247 void AliTRDrecoTask::Terminate(Option_t *)
248 {
249   //
250   // Terminate
251   //
252
253   if(fgDebugStream){ 
254     delete fgDebugStream;
255     fgDebugStream = NULL;
256   }
257   if(HasPostProcess()) PostProcess();
258 }
259
260 //________________________________________________________
261 void AliTRDrecoTask::Adjust(TF1 *f, TH1 * const h)
262 {
263 // Helper function to avoid duplication of code
264 // Make first guesses on the fit parameters
265
266   // find the intial parameters of the fit !! (thanks George)
267   Int_t nbinsy = Int_t(.5*h->GetNbinsX());
268   Double_t sum = 0.;
269   for(Int_t jbin=nbinsy-4; jbin<=nbinsy+4; jbin++) sum+=h->GetBinContent(jbin); sum/=9.;
270   f->SetParLimits(0, 0., 3.*sum);
271   f->SetParameter(0, .9*sum);
272
273   f->SetParLimits(1, -.2, .2);
274   f->SetParameter(1, -0.1);
275
276   f->SetParLimits(2, 0., 4.e-1);
277   f->SetParameter(2, 2.e-2);
278   if(f->GetNpar() <= 4) return;
279
280   f->SetParLimits(3, 0., sum);
281   f->SetParameter(3, .1*sum);
282
283   f->SetParLimits(4, -.3, .3);
284   f->SetParameter(4, 0.);
285
286   f->SetParLimits(5, 0., 1.e2);
287   f->SetParameter(5, 2.e-1);
288 }