]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/TRDcalib/AliTRDcalib.cxx
minor correction in documentation
[u/mrichter/AliRoot.git] / TRD / TRDcalib / AliTRDcalib.cxx
CommitLineData
da7d31b7 1#define AliTRDcalib_cxx
2// The class definition in Calib.h has been generated automatically
3// by the ROOT utility TTree::MakeSelector(). This class is derived
4// from the ROOT class TSelector. For more information on the TSelector
5// framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
6
7// The following methods are defined in this file:
8// Begin(): called everytime a loop on the tree starts,
9// a convenient place to create your histograms.
10// SlaveBegin(): called after Begin(), when on PROOF called only on the
11// slave servers.
12// Process(): called for each event, in this function you decide what
13// to read and fill your histograms.
14// SlaveTerminate: called at the end of the loop on the tree, when on PROOF
15// called only on the slave servers.
16// Terminate(): called at the end of the loop on the tree,
17// a convenient place to draw/fit your histograms.
18//
19// To use this file, try the following session on your Tree T:
20//
21// Root > T->Process("AliTRDcalib.cxx")
22// Root > T->Process("AliTRDcalib.cxx","some options")
23// Root > T->Process("AliTRDcalib.cxx+")
24//
25
26#include "AliTRDcalib.h"
27#include <TTree.h>
28#include <TObject.h>
29#include <TH2.h>
30#include <TStyle.h>
31#include <TH2I.h>
32#include <TProfile2D.h>
33#include <TCanvas.h>
34#include <TStyle.h>
35
36#include <AliESDEvent.h>
37#include <AliESDfriend.h>
38#include <AliESDtrack.h>
39#include <AliESDfriendTrack.h>
40
41#include <AliTRDgeometry.h>
42#include <AliTRDtrack.h>
43#include <AliTRDcluster.h>
44#include <AliCDBManager.h>
45#include <AliTRDCalibraFillHisto.h>
46#include <AliTRDCalibraVdriftLinearFit.h>
47
48
49AliTRDcalib::AliTRDcalib(TTree *) :
50 TSelector(),
51 fev(0),
52 fevf(0),
53 fo(0),
54 ft(0),
55 fc(0),
56 fesdTrack(0),
57 ffriendTrack(0),
58 fFileNo(0)
59 {
60 //G__SetCatchException(0);
61 }
62//_____________________________________________________________________
63void AliTRDcalib::Begin(TTree * /*tree*/)
64{
65 // The Begin() function is called at the start of the query.
66 // When running with PROOF Begin() is only called on the client.
67 // The tree argument is deprecated (on PROOF 0 is passed).
68
69 //TString option = GetOption();
70
71}
72//______________________________________________________________________
73void AliTRDcalib::SlaveBegin(TTree * tree)
74{
75 // The SlaveBegin() function is called after the Begin() function.
76 // When running with PROOF SlaveBegin() is called on each slave server.
77 // The tree argument is deprecated (on PROOF 0 is passed).
78
79 //printf("Slave Begin\n");
80
81 //TString option = GetOption();
82 if(tree) Init(tree);
83
84 fo = 0x0;
85 ft = 0x0;
86 fc = 0x0;
87 fesdTrack = 0x0;
88 ffriendTrack = 0x0;
89
90 AliCDBManager *cdbManager = AliCDBManager::Instance();
91 cdbManager->SetDefaultStorage("local://$ALICE_ROOT");
92 //cdbManager->SetSpecificStorage("TRD/Calib/FEE","local:///u/bailhach/aliroot/database30head/");
93 cdbManager->SetRun(0);
94
95 // instance calibration
96 fcalib = AliTRDCalibraFillHisto::Instance();
97 fcalib->SetNz(0,0);
98 fcalib->SetNrphi(0,0);
99 fcalib->SetNz(1,0);
100 fcalib->SetNrphi(1,0);
101 fcalib->SetMITracking();
102 fcalib->SetHisto2d();
103 fcalib->SetVector2d();
104 fcalib->SetLinearFitterOn();
105 fcalib->SetLinearFitterDebugOn();
106 fcalib->SetCH2dOn();
107 fcalib->SetPH2dOn();
108 fcalib->SetPRF2dOn();
109 fcalib->Init2Dhistostrack();
110 //fcalib->SetDebugLevel(1);
111 fcalib->SetNumberClusters(14);
112
113}
114//__________________________________________________________________________
115void AliTRDcalib::CleanESD(){
116 //
117 //printf("CleanESD\n");
118 if (fev!=0){
119 delete fev;
120 fev = 0;
121 }
122 if (fevf!=0){
123 delete fevf;
124 fevf =0;
125 }
126}
127//_________________________________________________________________________________
128Bool_t AliTRDcalib::Process(Long64_t entry)
129{
130 // The Process() function is called for each entry in the tree (or possibly
131 // keyed object in the case of PROOF) to be processed. The entry argument
132 // specifies which entry in the currently loaded tree is to be processed.
133 // It can be passed to either Calib::GetEntry() or TBranch::GetEntry()
134 // to read either all or the required parts of the data. When processing
135 // keyed objects with PROOF, the object is already loaded and is available
136 // via the fObject pointer.
137 //
138 // This function should contain the "body" of the analysis. It can contain
139 // simple or elaborate selection criteria, run algorithms on the data
140 // of the event and typically fill histograms.
141 //
142 // The processing can be stopped by calling Abort().
143 //
144 // Use fStatus to set the return value of TTree::Process().
145 //
146 // The return value is currently not used.
147 //printf("process\n");
148
149 if (!fChain) return kFALSE;
150 //printf("process1\n");
151 Int_t nBytes;
152 Int_t nTRDcls = 0;
153 nBytes = fChain->GetTree()->GetEntry(entry);
154 //printf("There are %d bytes for these event\n",nBytes);
155 if (!fev || (nBytes == 0)) {
156 return kFALSE;
157 }
158 if(fev->GetAliESDOld()) fev->CopyFromOldESD();
159 //printf("process2\n");
160 Int_t ntr = fev->GetNumberOfTracks();
161 //printf("Tracks new %d\n",ntr);
162
163 if (!fevf || (fevf->GetNumberOfTracks()!=ntr)) {
164 return kFALSE;
165 }
166
167 if(ntr>0){
168
169 fev->SetESDfriend(fevf);
170
171 //printf("Number of friends tracks %d\n",fevf->GetNumberOfTracks());
172
173
174 for(int itrk=0; itrk<fev->GetNumberOfTracks(); itrk++){
175
176 fesdTrack = fev->GetTrack(itrk);
177 if(!(nTRDcls = fesdTrack->GetTRDncls())) continue;
178 if(!(fesdTrack->GetFriendTrack())) continue;
179 ffriendTrack = new AliESDfriendTrack(*(fesdTrack->GetFriendTrack()));
180
181 Int_t icalib=0;
182 while((fo = (TObject *)(ffriendTrack->GetCalibObject(icalib++)))){
183 //printf("Name of calibObject %s\n",fo->IsA()->GetName());
184 if(strcmp(fo->IsA()->GetName(), "AliTRDtrack") != 0) continue;
185 //printf("\tfound %s @ 0x%x; calib object %d\n", fo->IsA()->GetName(), fo, icalib-1);
186 ft = (AliTRDtrack *)fo;
187
188 fcalib->UpdateHistograms(ft);
189 }
190 }
191 }
192
193 //CleanESD();
194 return kTRUE;
195}
196//______________________________________________________________________________________________
197void AliTRDcalib::SlaveTerminate()
198{
199 // The SlaveTerminate() function is called after all entries or objects
200 // have been processed. When running with PROOF SlaveTerminate() is called
201 // on each slave server.
202
203 if(!fOutput)
204 {
205 printf("ERROR: Output list not initialized\n");
206 return;
207 }
208
209 fCH2d = new TH2I(*(fcalib->GetCH2d()));
210 fPH2d = new TProfile2D(*(fcalib->GetPH2d()));
211 fPRF2d = new TProfile2D(*(fcalib->GetPRF2d()));
212
213 AliTRDCalibraVdriftLinearFit *ju = fcalib->GetVdriftLinearFit();
214 for(Int_t det = 0; det < 540; det++){
215 fVdriftLinear[det] = new TH2F(*(ju->GetLinearFitterHisto(det,kTRUE)));
216 }
217
218
219
220 fOutput->Add(fCH2d);
221 fOutput->Add(fPH2d);
222 fOutput->Add(fPRF2d);
223 for(Int_t det = 0; det < 540; det++){
224 fOutput->Add(fVdriftLinear[det]);
225 }
226
227 fcalib->Destroy();
228
229}
230//____________________________________________________________________________________
231void AliTRDcalib::Terminate()
232{
233 // The Terminate() function is the last function to be called during
234 // a query. It always runs on the client, it can be used to present
235 // the results graphically or save the results to file.
236
237
238 printf("InTerminate()\n");
239 if (!fOutput) return;
240
241 //fOutput->Print();
242
243 fCH2d = dynamic_cast<TH2I*>(fOutput->FindObject("CH2d"));
244
245 if (!fCH2d)
246 {
247 printf("Error: %s not returned\n","fCH2d");
248 return;
249 }
250
251 fPH2d = dynamic_cast<TProfile2D*>(fOutput->FindObject("PH2d"));
252
253 if (!fPH2d)
254 {
255 printf("Error: %s not returned\n","fPH2d");
256 return;
257 }
258
259 fPRF2d = dynamic_cast<TProfile2D*>(fOutput->FindObject("PRF2d"));
260
261 if (!fPRF2d)
262 {
263 printf("Error: %s not returned\n","fPRF2d");
264 return;
265 }
266
267
268 const char * Name = 0x0;
269 Name = "LFDV%dversion0";
270 TObjArray array(540);
271
272 for(Int_t det = 0; det < 540; det++){
273
274 TString Namehisto (Form(Name, det));
275 fVdriftLinear[det] = dynamic_cast<TH2F*>(fOutput->FindObject((const char *)Namehisto));
276 if (!fVdriftLinear[det])
277 {
278 printf("Error: %s not returned\n",(const char *)Namehisto);
279 return;
280 }
281 array.AddAt(fVdriftLinear[det],det);
282 }
283
284 AliTRDCalibraVdriftLinearFit vdriftlinearfit = AliTRDCalibraVdriftLinearFit(array);
285 vdriftlinearfit.FillPEArray();
286
287 gStyle->SetPalette(1);
288 gStyle->SetOptStat(1111);
289 gStyle->SetPadBorderMode(0);
290 gStyle->SetCanvasColor(10);
291 gStyle->SetPadLeftMargin(0.13);
292 gStyle->SetPadRightMargin(0.10);
293
294
295 printf("There are %d files analysed\n",fFileNo);
296
297
298 TCanvas *u = new TCanvas("u","",50,50,600,800);
299 u->Divide(3,1);
300 u->cd(1);
301 fCH2d->DrawCopy("lego");
302 u->cd(2);
303 fPH2d->DrawCopy("lego");
304 u->cd(3);
305 fPRF2d->DrawCopy("lego");
306
307 TObjArray *arrayE = vdriftlinearfit.GetEArray();
308
309 Double_t totalsum = 0.0;
310 for(Int_t k = 0; k < 540; k++){
311 TVectorD *h = (TVectorD *) arrayE->UncheckedAt(k);
312 if(h){
313 totalsum += (*h)[2];
314 }
315 }
316
317 TFile file("Output.root","recreate");
318 fOutput->Write();
319
320}
321//___________________________________________________________________________________________
322void AliTRDcalib::Init(TTree *tree)
323{
324 // The Init() function is called when the selector needs to initialize
325 // a new tree or chain. Typically here the branch addresses and branch
326 // pointers of the tree will be set.
327 // It is normaly not necessary to make changes to the generated
328 // code, but the routine can be extended by the user if needed.
329 // Init() will be called many times when running on PROOF
330 // (once per file to be processed).
331
332 // Set branch addresses and branch pointers
333 if (!tree) return;
334 fChain = tree;
335 fev = new AliESDEvent();
336 fChain->SetBranchStatus("ESDfriend*",1);
337 fev->ReadFromTree(fChain);
338 fevf = 0x0;
339 fevf = (AliESDfriend*)fev->FindListObject("AliESDfriend");
340 if(!fevf){
341 fChain->SetBranchAddress("ESDfriend.",&fevf);
342 }
343
344
345}
346//___________________________________________________________________________________________
347Bool_t AliTRDcalib::Notify()
348{
349 // The Notify() function is called when a new file is opened. This
350 // can be either for a new TTree in a TChain or when when a new TTree
351 // is started when using PROOF. It is normaly not necessary to make changes
352 // to the generated code, but the routine can be extended by the
353 // user if needed. The return value is currently not used.
354
355 ++fFileNo;
356 printf ("Processing file no %d\n",fFileNo);
357
358 return kTRUE;
359}