]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDtrackingResolution.cxx
- fix memory leaks
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDtrackingResolution.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercialf purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 /* $Id: AliTRDtrackingResolution.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  TRD tracking resolution                                               //
21 //
22 // The class performs resolution and residual studies 
23 // of the TRD tracks for the following quantities :
24 //   - spatial position (y, [z])
25 //   - angular (phi) tracklet
26 //   - momentum at the track level
27 // 
28 // The class has to be used for regular detector performance checks using the official macros:
29 //   - $ALICE_ROOT/TRD/qaRec/run.C
30 //   - $ALICE_ROOT/TRD/qaRec/makeResults.C
31 // 
32 // For stand alone usage please refer to the following example: 
33 // {  
34 //   gSystem->Load("libANALYSIS.so");
35 //   gSystem->Load("libTRDqaRec.so");
36 //   AliTRDtrackingResolution *res = new AliTRDtrackingResolution();
37 //   //res->SetMCdata();
38 //   //res->SetVerbose();
39 //   //res->SetVisual();
40 //   res->Load("TRD.TaskResolution.root");
41 //   if(!res->PostProcess()) return;
42 //   res->GetRefFigure(0);
43 // }  
44 //
45 //  Authors:                                                              //
46 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
47 //    Markus Fasel <M.Fasel@gsi.de>                                       //
48 //                                                                        //
49 ////////////////////////////////////////////////////////////////////////////
50
51 #include <cstring>
52
53 #include <TROOT.h>
54 #include <TSystem.h>
55 #include <TObjArray.h>
56 #include <TH2.h>
57 #include <TH1.h>
58 #include <TF1.h>
59 #include <TCanvas.h>
60 #include <TBox.h>
61 #include <TProfile.h>
62 #include <TGraphErrors.h>
63 #include <TMath.h>
64 #include "TTreeStream.h"
65 #include "TGeoManager.h"
66
67 #include "AliAnalysisManager.h"
68 #include "AliTrackReference.h"
69 #include "AliTrackPointArray.h"
70 #include "AliCDBManager.h"
71
72 #include "AliTRDSimParam.h"
73 #include "AliTRDgeometry.h"
74 #include "AliTRDpadPlane.h"
75 #include "AliTRDcluster.h"
76 #include "AliTRDseedV1.h"
77 #include "AliTRDtrackV1.h"
78 #include "AliTRDtrackerV1.h"
79 #include "AliTRDReconstructor.h"
80 #include "AliTRDrecoParam.h"
81
82 #include "AliTRDtrackInfo/AliTRDclusterInfo.h"
83 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
84 #include "AliTRDtrackingResolution.h"
85
86 ClassImp(AliTRDtrackingResolution)
87
88 //________________________________________________________
89 AliTRDtrackingResolution::AliTRDtrackingResolution()
90   :AliTRDrecoTask("Resolution", "Tracking Resolution")
91   ,fStatus(0)
92   ,fReconstructor(0x0)
93   ,fGeo(0x0)
94   ,fGraphS(0x0)
95   ,fGraphM(0x0)
96   ,fClResiduals(0x0)
97   ,fTrkltResiduals(0x0)
98   ,fTrkltPhiResiduals(0x0)
99   ,fClResolution(0x0)
100   ,fTrkltResolution(0x0)
101 {
102   fReconstructor = new AliTRDReconstructor();
103   fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
104   fGeo = new AliTRDgeometry();
105
106   InitFunctorList();
107
108   DefineOutput(1+kClusterResidual, TObjArray::Class());
109   DefineOutput(1+kTrackletYResidual, TObjArray::Class());
110   DefineOutput(1+kTrackletPhiResidual, TObjArray::Class());
111   DefineOutput(1+kClusterResolution, TObjArray::Class());
112   DefineOutput(1+kTrackletYResolution, TObjArray::Class());
113 }
114
115 //________________________________________________________
116 AliTRDtrackingResolution::~AliTRDtrackingResolution()
117 {
118   if(fGraphS){fGraphS->Delete(); delete fGraphS;}
119   if(fGraphM){fGraphM->Delete(); delete fGraphM;}
120   delete fGeo;
121   delete fReconstructor;
122   if(gGeoManager) delete gGeoManager;
123   if(fClResiduals){fClResiduals->Delete(); delete fClResiduals;}
124   if(fTrkltResiduals){fTrkltResiduals->Delete(); delete fTrkltResiduals;}
125   if(fTrkltPhiResiduals){fTrkltPhiResiduals->Delete(); delete fTrkltPhiResiduals;}
126   if(fClResolution){
127     fClResolution->Delete(); 
128     delete fClResolution;
129   }
130   if(fTrkltResolution){fTrkltResolution->Delete(); delete fTrkltResolution;}
131 }
132
133
134 //________________________________________________________
135 void AliTRDtrackingResolution::CreateOutputObjects()
136 {
137   // spatial resolution
138   OpenFile(0, "RECREATE");
139
140   fContainer = Histos();
141
142   fClResiduals = new TObjArray();
143   fClResiduals->SetOwner(kTRUE);
144   fTrkltResiduals = new TObjArray();
145   fTrkltResiduals->SetOwner(kTRUE);
146   fTrkltPhiResiduals = new TObjArray();
147   fTrkltPhiResiduals->SetOwner(kTRUE);
148   fClResolution = new TObjArray();
149   fClResolution->SetOwner(kTRUE);
150   fTrkltResolution = new TObjArray();
151   fTrkltResolution->SetOwner(kTRUE);
152 }
153
154 //________________________________________________________
155 void AliTRDtrackingResolution::Exec(Option_t *opt)
156 {
157   fClResiduals->Delete();
158   fTrkltResiduals->Delete();
159   fTrkltPhiResiduals->Delete();
160   fClResolution->Delete();
161   fTrkltResolution->Delete();
162
163   AliTRDrecoTask::Exec(opt);
164
165   PostData(1+kClusterResidual, fClResiduals);
166   PostData(1+kTrackletYResidual, fTrkltResiduals);
167   PostData(1+kTrackletPhiResidual, fTrkltPhiResiduals);
168   PostData(1+kClusterResolution, fClResolution);
169   PostData(1+kTrackletYResolution, fTrkltResolution);
170 }
171
172 //________________________________________________________
173 TH1* AliTRDtrackingResolution::PlotClusterResiduals(const AliTRDtrackV1 *track)
174 {
175   if(track) fTrack = track;
176   if(!fTrack){
177     AliWarning("No Track defined.");
178     return 0x0;
179   }
180   TH1 *h = 0x0;
181   if(!(h = ((TH2I*)fContainer->At(kClusterResidual)))){
182     AliWarning("No output histogram defined.");
183     return 0x0;
184   }
185
186   Float_t x0, y0, z0, dy, dydx, dzdx;
187   AliTRDseedV1 *fTracklet = 0x0;  
188   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
189     if(!(fTracklet = fTrack->GetTracklet(ily))) continue;
190     if(!fTracklet->IsOK()) continue;
191     x0 = fTracklet->GetX0();
192
193     // retrive the track angle with the chamber
194     y0   = fTracklet->GetYref(0);
195     z0   = fTracklet->GetZref(0);
196     dydx = fTracklet->GetYref(1);
197     dzdx = fTracklet->GetZref(1);
198     Float_t tilt = fTracklet->GetTilt();
199     AliTRDcluster *c = 0x0;
200     fTracklet->ResetClusterIter(kFALSE);
201     while((c = fTracklet->PrevCluster())){
202       Float_t xc = c->GetX();
203       Float_t yc = c->GetY();
204       Float_t zc = c->GetZ();
205       Float_t dx = x0 - xc; 
206       Float_t yt = y0 - dx*dydx;
207       Float_t zt = z0 - dx*dzdx; 
208       dy = yt - (yc - tilt*(zc-zt));
209
210       //dy = trklt.GetYat(c->GetX()) - c->GetY();
211       h->Fill(dydx, dy);
212   
213       if(fDebugLevel>=1){
214         // Get z-position with respect to anode wire
215         AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
216         Int_t istk = fGeo->GetStack(c->GetDetector());
217         AliTRDpadPlane *pp = fGeo->GetPadPlane(ily, istk);
218         Float_t row0 = pp->GetRow0();
219         Float_t d  =  row0 - zt + simParam->GetAnodeWireOffset();
220         d -= ((Int_t)(2 * d)) / 2.0;
221         if (d > 0.25) d  = 0.5 - d;
222
223 /*        AliTRDclusterInfo *clInfo = new AliTRDclusterInfo;
224         fClResiduals->Add(clInfo);
225         clInfo->SetCluster(c);
226         clInfo->SetGlobalPosition(yt, zt, dydx, dzdx);
227         clInfo->SetResolution(dy);
228         clInfo->SetAnisochronity(d);
229         clInfo->SetDriftLength(dx);
230         (*fDebugStream) << "ClusterResiduals"
231           <<"clInfo.=" << clInfo
232           << "\n";*/
233       }
234     }
235   }
236   return h;
237 }
238
239
240 //________________________________________________________
241 TH1* AliTRDtrackingResolution::PlotTrackletResiduals(const AliTRDtrackV1 *track)
242 {
243   if(track) fTrack = track;
244   if(!fTrack){
245     AliWarning("No Track defined.");
246     return 0x0;
247   }
248   TH1 *h = 0x0;
249   if(!(h = ((TH2I*)fContainer->At(kTrackletYResidual)))){
250     AliWarning("No output histogram defined.");
251     return 0x0;
252   }
253
254   AliTRDseedV1 *tracklet = 0x0;
255   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
256     if(!(tracklet = fTrack->GetTracklet(il))) continue;
257     if(!tracklet->IsOK()) continue;
258     h->Fill(tracklet->GetYref(1), tracklet->GetYref(0)-tracklet->GetYfit(0));
259   }
260   return h;
261
262   // refit the track
263 //   AliRieman fRim(fTrack->GetNumberOfClusters());
264 //   Float_t x[AliTRDgeometry::kNlayer] = {-1., -1., -1., -1., -1., -1.}, y[AliTRDgeometry::kNlayer], dydx[AliTRDgeometry::kNlayer];
265
266 //     AliTRDcluster *c = 0x0;
267 //     tracklet->ResetClusterIter(kFALSE);
268 //     while((c = tracklet->PrevCluster())){
269 //       Float_t xc = c->GetX();
270 //       Float_t yc = c->GetY();
271 //       Float_t zc = c->GetZ();
272 //       Float_t zt = tracklet->GetZref(0) - (tracklet->GetX0()-xc)*tracklet->GetZref(1); 
273 //       yc -= tracklet->GetTilt()*(zc-zt);
274 //       fRim.AddPoint(xc, yc, zc, 1, 10);
275 //     }
276 //     tracklet->Fit(kTRUE);
277 // 
278 //     x[il] = tracklet->GetX0();
279 //     y[il] = tracklet->GetYfit(0)-tracklet->GetYfit(1)*(tracklet->GetX0()-x[il]);
280 //     dydx[il] = tracklet->GetYref(1);
281
282
283 //   fRim.Update();
284
285 //   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
286 //     if(x[il] < 0.) continue;
287 //     Float_t dy = y[il]-fRim.GetYat(x[il])/*/sigma_track*/;
288 //     h->Fill(dydx[il], dy);
289 // 
290 //     if(fDebugLevel>=1){
291 //       Double_t sigmay = fRim.GetErrY(x[il]);
292 //       Float_t yt = fRim.GetYat(x[il]);
293 //       (*fDebugStream) << "TrkltResiduals"
294 //         << "layer="  << il
295 //         << "x="      << x[il]
296 //         << "y="      << y[il]
297 //         << "yt="     << yt
298 //         << "dydx="   << dydx[il]
299 //         << "dy="     << dy
300 //         << "sigmay=" << sigmay
301 //         << "\n";
302 //     }
303 //   }
304 }
305
306 //________________________________________________________
307 TH1* AliTRDtrackingResolution::PlotTrackletPhiResiduals(const AliTRDtrackV1 *track)
308 {
309   if(track) fTrack = track;
310   if(!fTrack){
311     AliWarning("No Track defined.");
312     return 0x0;
313   }
314   TH1 *h = 0x0;
315   if(!(h = ((TH2I*)fContainer->At(kTrackletPhiResidual)))){
316     AliWarning("No output histogram defined.");
317     return 0x0;
318   }
319
320   AliTRDseedV1 *tracklet = 0x0;
321   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
322     if(!(tracklet = fTrack->GetTracklet(il))) continue;
323     if(!tracklet->IsOK()) continue;
324     h->Fill(tracklet->GetYref(1), TMath::ATan(tracklet->GetYref(1))-TMath::ATan(tracklet->GetYfit(1)));
325   }
326   return h;
327
328 //   // refit the track
329 //   AliRieman fRim(fTrack->GetNumberOfClusters());
330 //   Float_t x[AliTRDgeometry::kNlayer] = {-1., -1., -1., -1., -1., -1.}, y[AliTRDgeometry::kNlayer], dydx[AliTRDgeometry::kNlayer];
331 //   Float_t dydx_ref=0, dydx_fit=0, phiref=0, phifit=0, phidiff=0;
332 //   AliTRDseedV1 *tracklet = 0x0;
333 //   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
334 //     if(!(tracklet = fTrack->GetTracklet(il))) continue;
335 //     if(!tracklet->IsOK()) continue;
336 //     AliTRDcluster *c = 0x0;
337 //     tracklet->ResetClusterIter(kFALSE);
338 //     while((c = tracklet->PrevCluster())){
339 //       Float_t xc = c->GetX();
340 //       Float_t yc = c->GetY();
341 //       Float_t zc = c->GetZ();
342 //       Float_t zt = tracklet->GetZref(0) - (tracklet->GetX0()-xc)*tracklet->GetZref(1); 
343 //       yc -= tracklet->GetTilt()*(zc-zt);
344 //       fRim.AddPoint(xc, yc, zc, 1, 10);
345 //     }
346 //     tracklet->Fit(kTRUE);
347 // 
348 //     x[il] = tracklet->GetX0();
349 //     y[il] = tracklet->GetYfit(0)-tracklet->GetYfit(1)*(tracklet->GetX0()-x[il]);
350 //     dydx[il] = tracklet->GetYref(1);
351 //   }
352 //   fRim.Update();
353 // 
354 //   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
355 //     if(x[il] < 0.) continue;
356 //     if(!(tracklet = fTrack->GetTracklet(il))) continue;
357 //     if(!tracklet->IsOK()) continue;
358 // 
359 //     dydx_ref = fRim.GetDYat(x[il]); 
360 //     dydx_fit = tracklet->GetYfit(1);
361 // 
362 //     phiref = TMath::ATan(dydx_ref);//*TMath::RadToDeg();
363 //     phifit = TMath::ATan(dydx_fit);//*TMath::RadToDeg();
364 //     
365 //     phidiff = phiref-phifit; /*/sigma_phi*/;
366 // 
367 //     h->Fill(dydx_ref, phidiff);
368 // 
369 // 
370 //     if(fDebugLevel>=1){
371 //       (*fDebugStream) << "TrkltPhiResiduals"
372 //         << "layer="  << il
373 //         << "dydx_fit="   << dydx_fit
374 //         << "dydx_ref="   << dydx_ref
375 //         << "phiref="     << phiref
376 //         << "phifit="     << phifit
377 //         << "phidiff="     << phidiff
378 //         << "\n";
379 //     }
380 //   }
381 }
382
383
384 //________________________________________________________
385 TH1* AliTRDtrackingResolution::PlotResolution(const AliTRDtrackV1 *track)
386 {
387   if(!HasMCdata()){ 
388     AliWarning("No MC defined. Results will not be available.");
389     return 0x0;
390   }
391   if(track) fTrack = track;
392   if(!fTrack){
393     AliWarning("No Track defined.");
394     return 0x0;
395   }
396   TH1 *h = 0x0;
397   if(!(h = ((TH2I*)fContainer->At(kClusterResolution)))){
398     AliWarning("No output histogram defined.");
399     return 0x0;
400   }
401   if(!(h = ((TH2I*)fContainer->At(kTrackletYResolution)))){
402     AliWarning("No output histogram defined.");
403     return 0x0;
404   }
405   if(!(h = ((TH2I*)fContainer->At(kTrackletZResolution)))){
406     AliWarning("No output histogram defined.");
407     return 0x0;
408   }
409   if(!(h = ((TH2I*)fContainer->At(kTrackletAngleResolution)))){
410     AliWarning("No output histogram defined.");
411     return 0x0;
412   }
413   //printf("called for %d tracklets ... \n", fTrack->GetNumberOfTracklets());
414   UChar_t s;
415   Int_t pdg = fMC->GetPDG(), det=-1;
416   Int_t label = fMC->GetLabel();
417   Float_t x0, y0, z0, dx, dy, dydx, dzdx;
418   AliTRDseedV1 *fTracklet = 0x0;  
419   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
420     if(!(fTracklet = fTrack->GetTracklet(ily))) continue;
421     if(!fTracklet->IsOK()) continue;
422     //printf("process layer[%d]\n", ily);
423     // retrive the track position and direction within the chamber
424     det = fTracklet->GetDetector();
425     x0  = fTracklet->GetX0();
426     if(!fMC->GetDirections(x0, y0, z0, dydx, dzdx, s)) continue; 
427
428     // recalculate tracklet based on the MC info
429     AliTRDseedV1 tt(*fTracklet);
430     tt.SetZref(0, z0);
431     tt.SetZref(1, dzdx); 
432     if(!tt.Fit(kTRUE)) continue;
433     //tt.Update();
434
435     dx = 0.;//x0 - tt.GetXref();
436     Float_t yt = y0 - dx*dydx;
437     Float_t yf = tt.GetYfit(0) - (fTracklet->GetX0() - x0)*tt.GetYfit(1);
438     dy = yf-yt;
439     Float_t dphi   = TMath::ATan(tt.GetYfit(1)) - TMath::ATan(dydx);
440     Float_t dz = 100.;
441     Bool_t  cross = fTracklet->GetNChange();
442     if(cross){
443       Double_t *xyz = tt.GetCrossXYZ();
444       dz = xyz[2] - (z0 - (x0 - xyz[0])*dzdx) ;
445       ((TH2I*)fContainer->At(kTrackletZResolution))->Fill(dzdx, dz);
446     }
447   
448     // Fill Histograms
449     ((TH2I*)fContainer->At(kTrackletYResolution))->Fill(dydx, dy);
450     ((TH2I*)fContainer->At(kTrackletAngleResolution))->Fill(dydx, dphi*TMath::RadToDeg());
451
452     // Fill Debug stream
453     if(fDebugLevel>=1){
454       Float_t p = fMC->GetTrackRef() ? fMC->GetTrackRef()->P() : -1.;
455       (*fDebugStream) << "TrkltResolution"
456         << "det="                 << det
457         << "pdg="     << pdg
458         << "p="       << p
459         << "ymc="     << y0
460         << "zmc="     << z0
461         << "dydx="    << dydx
462         << "dzdx="    << dzdx
463         << "cross="   << cross
464         << "dy="                  << dy
465         << "dz="                  << dz
466         << "dphi="              << dphi
467         << "\n";
468     }
469
470     Int_t istk = AliTRDgeometry::GetStack(det); 
471     AliTRDpadPlane *pp = fGeo->GetPadPlane(ily, istk);
472     Float_t zr0 = pp->GetRow0() + AliTRDSimParam::Instance()->GetAnodeWireOffset();
473     Float_t tilt = fTracklet->GetTilt();
474
475     AliTRDcluster *c = 0x0;
476     fTracklet->ResetClusterIter(kFALSE);
477     while((c = fTracklet->PrevCluster())){
478       Float_t  q = TMath::Abs(c->GetQ());
479       Float_t xc = c->GetX();
480       Float_t yc = c->GetY();
481       Float_t zc = c->GetZ();
482       dx = x0 - xc; 
483       Float_t yt = y0 - dx*dydx;
484       Float_t zt = z0 - dx*dzdx; 
485       dy = yt - (yc - tilt*(zc-zt));
486
487       // Fill Histograms
488       if(q>20. && q<250.) ((TH2I*)fContainer->At(kClusterResolution))->Fill(dydx, dy);
489       
490       // Fill Debug Tree
491       if(fDebugLevel>=1){
492         Float_t d = zr0 - zt;
493         d -= ((Int_t)(2 * d)) / 2.0;
494         if (d > 0.25) d  = 0.5 - d;
495
496         AliTRDclusterInfo *clInfo = new AliTRDclusterInfo;
497         fClResolution->Add(clInfo);
498         clInfo->SetCluster(c);
499         clInfo->SetMC(pdg, label);
500         clInfo->SetGlobalPosition(yt, zt, dydx, dzdx);
501         clInfo->SetResolution(dy);
502         clInfo->SetAnisochronity(d);
503         clInfo->SetDriftLength(dx);
504         clInfo->SetTilt(tilt);
505         //clInfo->Print();
506         (*fDebugStream) << "ClusterResolution"
507           <<"clInfo.=" << clInfo
508           << "\n";
509       }
510     }
511   }
512   return h;
513 }
514
515
516 //________________________________________________________
517 void AliTRDtrackingResolution::GetRefFigure(Int_t ifig)
518 {
519   TBox *b = 0x0;
520   TAxis *ax = 0x0;
521   TGraphErrors *g = 0x0;
522   switch(ifig){
523   case kClusterResidual:
524     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
525     g->Draw("apl");
526     ax = g->GetHistogram()->GetYaxis();
527     ax->SetRangeUser(-.5, 2.5);
528     ax->SetTitle("Clusters Y Residuals #sigma/#mu [mm]");
529     ax = g->GetHistogram()->GetXaxis();
530     ax->SetTitle("tg(#phi)");
531     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
532     g->Draw("pl");
533     b = new TBox(-.15, -.5, .15, 1.);
534     b->SetFillStyle(3002);b->SetFillColor(kGreen);
535     b->SetLineColor(0); b->Draw();
536     return;
537   case kTrackletYResidual:
538     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
539     g->Draw("apl");
540     ax = g->GetHistogram()->GetYaxis();
541     ax->SetRangeUser(-.5, 3.);
542     ax->SetTitle("Tracklet Y Residuals #sigma/#mu [mm]");
543     ax = g->GetHistogram()->GetXaxis();
544     ax->SetTitle("tg(#phi)");
545     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
546     g->Draw("pl");
547     b = new TBox(-.15, -.5, .15, 1.);
548     b->SetFillStyle(3002);b->SetFillColor(kGreen);
549     b->SetLineColor(0); b->Draw();
550     return;
551   case kTrackletPhiResidual:
552     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
553     g->Draw("apl");
554     ax = g->GetHistogram()->GetYaxis();
555     ax->SetRangeUser(-.5, 2.);
556     ax->SetTitle("Tracklet Phi Residuals #sigma/#mu [rad]");
557     ax = g->GetHistogram()->GetXaxis();
558     ax->SetTitle("tg(#phi)");
559     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
560     g->Draw("pl");
561     b = new TBox(-.15, -.5, .15, 1.);
562     b->SetFillStyle(3002);b->SetFillColor(kGreen);
563     b->SetLineColor(0); b->Draw();
564     return;
565   case kClusterResolution:
566     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
567     ax = g->GetHistogram()->GetYaxis();
568     ax->SetRangeUser(-.5, 1.);
569     ax->SetTitle("Cluster Y Resolution #sigma/#mu [mm]");
570     ax = g->GetHistogram()->GetXaxis();
571     ax->SetTitle("tg(#phi)");
572     g->Draw("apl");
573     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
574     g->Draw("pl");
575     b = new TBox(-.15, -.5, .15, 1.);
576     b->SetFillStyle(3002);b->SetFillColor(kGreen);
577     b->SetLineColor(0); b->Draw();
578     return;
579   case kTrackletYResolution:
580     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
581     ax = g->GetHistogram()->GetYaxis();
582     ax->SetRangeUser(-.5, 1.);
583     ax->SetTitle("Tracklet Y Resolution #sigma/#mu [mm]");
584     ax = g->GetHistogram()->GetXaxis();
585     ax->SetTitle("tg(#phi)");
586     g->Draw("apl");
587     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
588     g->Draw("pl");
589     return;
590   case kTrackletZResolution:
591     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
592     ax = g->GetHistogram()->GetYaxis();
593     ax->SetRangeUser(-.5, 1.);
594     ax->SetTitle("Tracklet Z Resolution #sigma/#mu [mm]");
595     ax = g->GetHistogram()->GetXaxis();
596     ax->SetTitle("tg(#theta)");
597     g->Draw("apl");
598     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
599     g->Draw("pl");
600     return;
601   case kTrackletAngleResolution:
602     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
603     ax = g->GetHistogram()->GetYaxis();
604     ax->SetRangeUser(-.05, .2);
605     ax->SetTitle("Tracklet Angular Resolution #sigma/#mu [deg]");
606     ax = g->GetHistogram()->GetXaxis();
607     ax->SetTitle("tg(#phi)");
608     g->Draw("apl");
609     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
610     g->Draw("pl");
611     return;
612   default:
613     AliInfo(Form("Reference plot [%d] not implemented yet", ifig));
614     return;
615   }
616   AliInfo(Form("Reference plot [%d] missing result", ifig));
617 }
618
619
620 //________________________________________________________
621 Bool_t AliTRDtrackingResolution::PostProcess()
622 {
623   //fContainer = dynamic_cast<TObjArray*>(GetOutputData(0));
624   if (!fContainer) {
625     Printf("ERROR: list not available");
626     return kFALSE;
627   }
628   fNRefFigures = fContainer->GetEntriesFast();
629   if(!fGraphS){ 
630     fGraphS = new TObjArray(fNRefFigures);
631     fGraphS->SetOwner();
632   }
633   if(!fGraphM){ 
634     fGraphM = new TObjArray(fNRefFigures);
635     fGraphM->SetOwner();
636   }
637
638   TH2I *h2 = 0x0;
639   TH1D *h = 0x0;
640   TGraphErrors *gm = 0x0, *gs = 0x0;
641
642   // define models
643   TF1 f("f1", "gaus", -.5, .5);  
644
645   TF1 fb("fb", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]", -.5, .5);
646
647   TF1 fc("fc", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]*exp(-0.5*((x-[4])/[5])**2)", -.5, .5);
648
649   TCanvas *c = 0x0;
650   if(IsVisual()) c = new TCanvas("c", Form("%s Visual", GetName()), 500, 500);
651   char opt[5];
652   sprintf(opt, "%c%c", IsVerbose() ? ' ' : 'Q', IsVisual() ? ' ': 'N');
653
654
655   //PROCESS RESIDUAL DISTRIBUTIONS
656
657   // Clusters residuals
658   h2 = (TH2I *)(fContainer->At(kClusterResidual));
659   gm = new TGraphErrors();
660   gm->SetLineColor(kBlue);
661   gm->SetMarkerStyle(7);
662   gm->SetMarkerColor(kBlue);
663   gm->SetNameTitle("clm", "");
664   fGraphM->AddAt(gm, kClusterResidual);
665   gs = new TGraphErrors();
666   gs->SetLineColor(kRed);
667   gs->SetMarkerStyle(23);
668   gs->SetMarkerColor(kRed);
669   gs->SetNameTitle("cls", "");
670   fGraphS->AddAt(gs, kClusterResidual);
671   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
672     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
673     h = h2->ProjectionY("py", ibin, ibin);
674     if(h->GetEntries()<100) continue;
675     AdjustF1(h, &f);
676
677     if(IsVisual()){c->cd(); c->SetLogy();}
678     h->Fit(&f, opt, "", -0.5, 0.5);
679     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
680     
681     Int_t ip = gm->GetN();
682     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
683     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
684     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
685     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
686   }
687
688   // Tracklet y residuals
689   h2 = (TH2I *)(fContainer->At(kTrackletYResidual));
690   gm = new TGraphErrors();
691   gm->SetLineColor(kBlue);
692   gm->SetMarkerStyle(7);
693   gm->SetMarkerColor(kBlue);
694   gm->SetNameTitle("tktm", "");
695   fGraphM->AddAt(gm, kTrackletYResidual);
696   gs = new TGraphErrors();
697   gs->SetLineColor(kRed);
698   gs->SetMarkerStyle(23);
699   gs->SetMarkerColor(kRed);
700   gs->SetNameTitle("tkts", "");
701   fGraphS->AddAt(gs, kTrackletYResidual);
702   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
703     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
704     h = h2->ProjectionY("py", ibin, ibin);
705     if(h->GetEntries()<100) continue;
706     AdjustF1(h, &f);
707
708     if(IsVisual()){c->cd(); c->SetLogy();}
709     h->Fit(&f, opt, "", -0.5, 0.5);
710     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
711     
712     Int_t ip = gm->GetN();
713     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
714     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
715     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
716     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
717   }
718
719   // Tracklet phi residuals
720   h2 = (TH2I *)(fContainer->At(kTrackletPhiResidual));
721   gm = new TGraphErrors();
722   gm->SetLineColor(kBlue);
723   gm->SetMarkerStyle(7);
724   gm->SetMarkerColor(kBlue);
725   gm->SetNameTitle("tktphim", "");
726   fGraphM->AddAt(gm, kTrackletPhiResidual);
727   gs = new TGraphErrors();
728   gs->SetLineColor(kRed);
729   gs->SetMarkerStyle(23);
730   gs->SetMarkerColor(kRed);
731   gs->SetNameTitle("tktphis", "");
732   fGraphS->AddAt(gs, kTrackletPhiResidual);
733   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
734     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
735     h = h2->ProjectionY("py", ibin, ibin);
736     if(h->GetEntries()<100) continue;
737     AdjustF1(h, &f);
738
739     if(IsVisual()){c->cd(); c->SetLogy();}
740     h->Fit(&f, opt, "", -0.5, 0.5);
741     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
742     
743     Int_t ip = gm->GetN();
744     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
745     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
746     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
747     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
748   }
749
750
751   //PROCESS RESOLUTION DISTRIBUTIONS
752
753   if(HasMCdata()){
754     // cluster y resolution
755     h2 = (TH2I*)fContainer->At(kClusterResolution);
756     gm = new TGraphErrors();
757     gm->SetLineColor(kBlue);
758     gm->SetMarkerStyle(7);
759     gm->SetMarkerColor(kBlue);
760     gm->SetNameTitle("clym", "");
761     fGraphM->AddAt(gm, kClusterResolution);
762     gs = new TGraphErrors();
763     gs->SetLineColor(kRed);
764     gs->SetMarkerStyle(23);
765     gs->SetMarkerColor(kRed);
766     gs->SetNameTitle("clys", "");
767     fGraphS->AddAt(gs, kClusterResolution);
768     for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
769       h = h2->ProjectionY("py", iphi, iphi);
770       if(h->GetEntries()<100) continue;
771       AdjustF1(h, &f);
772
773       if(IsVisual()){c->cd(); c->SetLogy();}
774       h->Fit(&f, opt, "", -0.5, 0.5);
775       if(IsVerbose()){
776         printf("phi[%d] mean[%e] sigma[%e]\n\n", iphi, 10.*f.GetParameter(1), 10.*f.GetParameter(2));
777       }
778       if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
779
780       Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
781       Int_t ip = gm->GetN();
782       gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
783       gm->SetPointError(ip, 0., 10.*f.GetParError(1));
784       gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
785       gs->SetPointError(ip, 0., 10.*f.GetParError(2));
786     }
787   
788     // tracklet y resolution
789     h2 = (TH2I*)fContainer->At(kTrackletYResolution);
790     gm = new TGraphErrors();
791     gm->SetLineColor(kBlue);
792     gm->SetMarkerStyle(7);
793     gm->SetMarkerColor(kBlue);
794     gm->SetNameTitle("trkltym", "");
795     fGraphM->AddAt(gm, kTrackletYResolution);
796     gs = new TGraphErrors();
797     gs->SetLineColor(kRed);
798     gs->SetMarkerStyle(23);
799     gs->SetMarkerColor(kRed);
800     gs->SetNameTitle("trkltys", "");
801     fGraphS->AddAt(gs, kTrackletYResolution);
802     for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
803       h = h2->ProjectionY("py", iphi, iphi);
804       if(h->GetEntries()<100) continue;
805       AdjustF1(h, &f);
806
807       if(IsVisual()){c->cd(); c->SetLogy();}
808       h->Fit(&f, opt, "", -0.5, 0.5);
809       if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
810
811       Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
812       Int_t ip = gm->GetN();
813       gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
814       gm->SetPointError(ip, 0., 10.*f.GetParError(1));
815       gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
816       gs->SetPointError(ip, 0., 10.*f.GetParError(2));
817     }
818   
819     // tracklet z resolution
820     h2 = (TH2I*)fContainer->At(kTrackletZResolution);
821     gm = new TGraphErrors();
822     gm->SetLineColor(kBlue);
823     gm->SetMarkerStyle(7);
824     gm->SetMarkerColor(kBlue);
825     gm->SetNameTitle("trkltzm", "");
826     fGraphM->AddAt(gm, kTrackletZResolution);
827     gs = new TGraphErrors();
828     gs->SetLineColor(kRed);
829     gs->SetMarkerStyle(23);
830     gs->SetMarkerColor(kRed);
831     gs->SetNameTitle("trkltzs", "");
832     fGraphS->AddAt(gs, kTrackletZResolution);
833     for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
834       h = h2->ProjectionY("py", iphi, iphi);
835       if(h->GetEntries()<100) continue;
836       AdjustF1(h, &fb);
837
838       if(IsVisual()){c->cd(); c->SetLogy();}
839       h->Fit(&fb, opt, "", -0.5, 0.5);
840       if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
841
842       Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
843       Int_t ip = gm->GetN();
844       gm->SetPoint(ip, phi, 10.*fb.GetParameter(1));
845       gm->SetPointError(ip, 0., 10.*fb.GetParError(1));
846       gs->SetPoint(ip, phi, 10.*fb.GetParameter(2));
847       gs->SetPointError(ip, 0., 10.*fb.GetParError(2));
848     }
849   
850     //tracklet phi resolution
851     h2 = (TH2I*)fContainer->At(kTrackletAngleResolution);
852     gm = new TGraphErrors();
853     gm->SetLineColor(kBlue);
854     gm->SetMarkerStyle(7);
855     gm->SetMarkerColor(kBlue);
856     gm->SetNameTitle("trkltym", "");
857     fGraphM->AddAt(gm, kTrackletAngleResolution);
858     gs = new TGraphErrors();
859     gs->SetLineColor(kRed);
860     gs->SetMarkerStyle(23);
861     gs->SetMarkerColor(kRed);
862     gs->SetNameTitle("trkltys", "");
863     fGraphS->AddAt(gs, kTrackletAngleResolution);
864     for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
865       h = h2->ProjectionY("py", iphi, iphi);
866       if(h->GetEntries()<100) continue;
867
868       if(IsVisual()){c->cd(); c->SetLogy();}
869       h->Fit(&f, opt, "", -0.5, 0.5);
870       if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
871
872       Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
873       Int_t ip = gm->GetN();
874       gm->SetPoint(ip, phi, f.GetParameter(1));
875       gm->SetPointError(ip, 0., f.GetParError(1));
876       gs->SetPoint(ip, phi, f.GetParameter(2));
877       gs->SetPointError(ip, 0., f.GetParError(2));
878     }
879   }
880   if(c) delete c;
881
882   return kTRUE;
883 }
884
885
886 //________________________________________________________
887 void AliTRDtrackingResolution::Terminate(Option_t *)
888 {
889   if(fDebugStream){ 
890     delete fDebugStream;
891     fDebugStream = 0x0;
892     fDebugLevel = 0;
893   }
894   if(HasPostProcess()) PostProcess();
895 }
896
897 //________________________________________________________
898 void AliTRDtrackingResolution::AdjustF1(TH1 *h, TF1 *f)
899 {
900 // Helper function to avoid duplication of code
901 // Make first guesses on the fit parameters
902
903   // find the intial parameters of the fit !! (thanks George)
904   Int_t nbinsy = Int_t(.5*h->GetNbinsX());
905   Double_t sum = 0.;
906   for(Int_t jbin=nbinsy-4; jbin<=nbinsy+4; jbin++) sum+=h->GetBinContent(jbin); sum/=9.;
907   f->SetParLimits(0, 0., 3.*sum);
908   f->SetParameter(0, .9*sum);
909
910   f->SetParLimits(1, -.2, .2);
911   f->SetParameter(1, -0.1);
912
913   f->SetParLimits(2, 0., 4.e-1);
914   f->SetParameter(2, 2.e-2);
915   if(f->GetNpar() <= 4) return;
916
917   f->SetParLimits(3, 0., sum);
918   f->SetParameter(3, .1*sum);
919
920   f->SetParLimits(4, -.3, .3);
921   f->SetParameter(4, 0.);
922
923   f->SetParLimits(5, 0., 1.e2);
924   f->SetParameter(5, 2.e-1);
925 }
926
927 //________________________________________________________
928 TObjArray* AliTRDtrackingResolution::Histos()
929 {
930   if(fContainer) return fContainer;
931
932   fContainer  = new TObjArray(7);
933   //fContainer->SetOwner(kTRUE);
934
935   TH1 *h = 0x0;
936   // cluster to tracklet residuals [2]
937   if(!(h = (TH2I*)gROOT->FindObject("fYCl"))){
938     h = new TH2I("fYCl", "Clusters Residuals", 21, -.33, .33, 100, -.5, .5);
939     h->GetXaxis()->SetTitle("tg(#phi)");
940     h->GetYaxis()->SetTitle("#Delta y [cm]");
941     h->GetZaxis()->SetTitle("entries");
942   } else h->Reset();
943   fContainer->AddAt(h, kClusterResidual);
944
945   // tracklet to track residuals [2]
946   if(!(h = (TH2I*)gROOT->FindObject("hTrkltYRez"))){ 
947     h = new TH2I("hTrkltYRez", "Tracklets", 21, -.33, .33, 100, -.5, .5);
948     h->GetXaxis()->SetTitle("tg(#phi)");
949     h->GetYaxis()->SetTitle("#Delta y [cm]");
950     h->GetZaxis()->SetTitle("entries");
951   } else h->Reset();
952   fContainer->AddAt(h, kTrackletYResidual);
953
954   // tracklet to track residuals angular [2]
955   if(!(h = (TH2I*)gROOT->FindObject("hTrkltPhiRez"))){ 
956     h = new TH2I("hTrkltPhiRez", "Tracklets", 21, -.33, .33, 100, -.5, .5);
957     h->GetXaxis()->SetTitle("tg(#phi)");
958     h->GetYaxis()->SetTitle("#Delta phi [#circ]");
959     h->GetZaxis()->SetTitle("entries");
960   } else h->Reset();
961   fContainer->AddAt(h, kTrackletPhiResidual);
962
963
964   // Resolution histos
965   if(HasMCdata()){
966     // cluster y resolution [0]
967     if(!(h = (TH2I*)gROOT->FindObject("fYClMC"))){
968       h = new TH2I("fYClMC", "Cluster Resolution", 31, -.48, .48, 100, -.5, .5);
969       h->GetXaxis()->SetTitle("tg(#phi)");
970       h->GetYaxis()->SetTitle("#Delta y [cm]");
971       h->GetZaxis()->SetTitle("entries");
972     } else h->Reset();
973     fContainer->AddAt(h, kClusterResolution);
974
975     // tracklet y resolution [0]
976     if(!(h = (TH2I*)gROOT->FindObject("fYTrkltMC"))){
977       h = new TH2I("fYTrkltMC", "Tracklet Resolution (Y)", 31, -.48, .48, 100, -.5, .5);
978       h->GetXaxis()->SetTitle("tg(#phi)");
979       h->GetYaxis()->SetTitle("#Delta y [cm]");
980       h->GetZaxis()->SetTitle("entries");
981     } else h->Reset();
982     fContainer->AddAt(h, kTrackletYResolution);
983
984     // tracklet y resolution [0]
985     if(!(h = (TH2I*)gROOT->FindObject("fZTrkltMC"))){
986       h = new TH2I("fZTrkltMC", "Tracklet Resolution (Z)", 31, -.48, .48, 100, -.5, .5);
987       h->GetXaxis()->SetTitle("tg(#theta)");
988       h->GetYaxis()->SetTitle("#Delta z [cm]");
989       h->GetZaxis()->SetTitle("entries");
990     } else h->Reset();
991     fContainer->AddAt(h, kTrackletZResolution);
992
993     // tracklet angular resolution [1]
994     if(!(h = (TH2I*)gROOT->FindObject("fPhiTrkltMC"))){
995       h = new TH2I("fPhiTrkltMC", "Tracklet Resolution (Angular)", 31, -.48, .48, 100, -10., 10.);
996       h->GetXaxis()->SetTitle("tg(#phi)");
997       h->GetYaxis()->SetTitle("#Delta #phi [deg]");
998       h->GetZaxis()->SetTitle("entries");
999     } else h->Reset();
1000     fContainer->AddAt(h, kTrackletAngleResolution);
1001   }
1002   return fContainer;
1003 }
1004
1005
1006 //________________________________________________________
1007 void AliTRDtrackingResolution::SetRecoParam(AliTRDrecoParam *r)
1008 {
1009
1010   fReconstructor->SetRecoParam(r);
1011 }