]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDrecoTask.cxx
fix ordering of tasks
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDrecoTask.cxx
1 #include "TClass.h"
2 #include "TMethod.h"
3 #include "TMethodCall.h"
4 #include "TMethodArg.h"
5 #include "TFile.h"
6 #include "TList.h"
7 #include "TH1.h"
8 #include "TF1.h"
9 #include "TObjArray.h"
10 #include "TDirectory.h"
11 #include "TTreeStream.h"
12
13 #include "AliLog.h"
14 #include "AliAnalysisTask.h"
15
16 #include "AliTRDrecoTask.h"
17
18 ClassImp(AliTRDrecoTask)
19
20 //_______________________________________________________
21 AliTRDrecoTask::AliTRDrecoTask(const char *name, const char *title)
22   : AliAnalysisTask(name, title)
23   ,fNRefFigures(0)
24   ,fDebugLevel(0)
25   ,fPlotFuncList(0x0)
26   ,fContainer(0x0)
27   ,fTracks(0x0)
28   ,fTrack(0x0)
29   ,fMC(0x0)
30   ,fESD(0x0)
31   ,fDebugStream(0x0)
32 {
33   DefineInput(0, TObjArray::Class());
34   DefineOutput(0, TObjArray::Class());
35 }
36
37 //_______________________________________________________
38 AliTRDrecoTask::~AliTRDrecoTask() 
39 {
40   printf(" %s (%s)\n", GetName(), GetTitle());
41   if(fDebugStream){ 
42     delete fDebugStream;
43     fDebugStream = 0x0;
44   }
45
46   if(fPlotFuncList){
47     fPlotFuncList->Delete();
48     delete fPlotFuncList;
49     fPlotFuncList = 0x0;
50   }
51   
52   if(fContainer){
53     if(fContainer->IsOwner()) fContainer->Delete();
54     delete fContainer;
55     fContainer = 0x0;
56   }
57 }
58
59 //_______________________________________________________
60 void AliTRDrecoTask::ConnectInputData(Option_t *)
61 {
62   //
63   // Connect input data
64   //
65
66   fTracks = dynamic_cast<TObjArray *>(GetInputData(0));
67 }
68
69 //_______________________________________________________
70 void AliTRDrecoTask::Exec(Option_t *)
71 {
72   if(!fPlotFuncList){
73     AliWarning("No functor list defined for the reference plots");
74     return;
75   }
76   if(!fTracks) return;
77   if(!fTracks->GetEntriesFast()) return;
78   
79   AliTRDtrackInfo *trackInfo = 0x0;
80   TIter plotIter(fPlotFuncList);
81   TObjArrayIter trackIter(fTracks);
82   while((trackInfo = dynamic_cast<AliTRDtrackInfo*>(trackIter()))){
83     fTrack = trackInfo->GetTrack();
84     fMC    = trackInfo->GetMCinfo();
85     fESD   = trackInfo->GetESDinfo();
86
87     TMethodCall *plot = 0x0;
88     plotIter.Reset();
89     while((plot=dynamic_cast<TMethodCall*>(plotIter()))){
90       plot->Execute(this);
91     }
92   }
93   PostData(0, fContainer);
94 }
95
96 //_______________________________________________________
97 Bool_t AliTRDrecoTask::GetRefFigure(Int_t /*ifig*/)
98 {
99   AliWarning("Retrieving reference figures not implemented.");
100   return kFALSE;
101 }
102
103 //_______________________________________________________
104 void AliTRDrecoTask::InitFunctorList()
105 {
106   TClass *c = this->IsA();
107
108   TMethod *m = 0x0;
109   TIter methIter(c->GetListOfMethods());
110   while((m=dynamic_cast<TMethod*>(methIter()))){
111     TString name(m->GetName());
112     if(!name.BeginsWith("Plot")) continue;
113     if(!fPlotFuncList) fPlotFuncList = new TList();
114     fPlotFuncList->AddLast(new TMethodCall(c, (const char*)name, ""));
115   }
116 }
117
118 //_______________________________________________________
119 Bool_t AliTRDrecoTask::Load(const Char_t *filename)
120 {
121   if(!TFile::Open(filename)){
122     AliWarning(Form("Couldn't open file %s.", filename));
123     return kFALSE;
124   }
125   TObjArray *o = 0x0;
126   if(!(o = (TObjArray*)gFile->Get(GetName()))){
127     AliWarning("Missing histogram container.");
128     return kFALSE;
129   }
130   fContainer = (TObjArray*)o->Clone(GetName());
131   gFile->Close();
132   return kTRUE;
133 }
134
135 //________________________________________________________
136 Bool_t AliTRDrecoTask::Save(TObjArray *results){
137   //
138   // Store the output graphs in a ROOT file
139   // Input TObject array will not be written as Key to the file,
140   // only content itself
141   //
142
143   TDirectory *cwd = gDirectory;
144   if(!TFile::Open(Form("TRD.Result%s.root", GetName()), "RECREATE")) return kFALSE;
145
146   TIterator *iter = results->MakeIterator();
147   TObject *inObject = 0x0, *outObject = 0x0;
148   while((inObject = iter->Next())){
149     outObject = inObject->Clone();
150     outObject->Write(0x0, TObject::kSingleKey);
151   }
152   delete iter;
153   gFile->Close(); delete gFile;
154   cwd->cd(); 
155   return kTRUE;
156 }
157
158 //_______________________________________________________
159 Bool_t AliTRDrecoTask::PostProcess()
160 {
161   AliWarning("Post processing of reference histograms not implemented.");
162   return kFALSE;
163 }
164
165 //_______________________________________________________
166 void AliTRDrecoTask::SetDebugLevel(Int_t level)
167 {
168   fDebugLevel = level;
169   if(fDebugLevel>=1){
170     TDirectory *savedir = gDirectory;
171     fDebugStream = new TTreeSRedirector(Form("TRD.Debug%s.root", GetName()));
172     savedir->cd();
173   }
174 }
175
176 //________________________________________________________
177 void AliTRDrecoTask::Adjust(TF1 *f, TH1 *h)
178 {
179 // Helper function to avoid duplication of code
180 // Make first guesses on the fit parameters
181
182   // find the intial parameters of the fit !! (thanks George)
183   Int_t nbinsy = Int_t(.5*h->GetNbinsX());
184   Double_t sum = 0.;
185   for(Int_t jbin=nbinsy-4; jbin<=nbinsy+4; jbin++) sum+=h->GetBinContent(jbin); sum/=9.;
186   f->SetParLimits(0, 0., 3.*sum);
187   f->SetParameter(0, .9*sum);
188
189   f->SetParLimits(1, -.2, .2);
190   f->SetParameter(1, -0.1);
191
192   f->SetParLimits(2, 0., 4.e-1);
193   f->SetParameter(2, 2.e-2);
194   if(f->GetNpar() <= 4) return;
195
196   f->SetParLimits(3, 0., sum);
197   f->SetParameter(3, .1*sum);
198
199   f->SetParLimits(4, -.3, .3);
200   f->SetParameter(4, 0.);
201
202   f->SetParLimits(5, 0., 1.e2);
203   f->SetParameter(5, 2.e-1);
204 }