1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDtrackingResolution.cxx 27496 2008-07-22 08:35:45Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // TRD tracking resolution //
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
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
32 // For stand alone usage please refer to the following example:
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);
46 // Alexandru Bercuci <A.Bercuci@gsi.de> //
47 // Markus Fasel <M.Fasel@gsi.de> //
49 ////////////////////////////////////////////////////////////////////////////
55 #include <TObjArray.h>
62 #include <TGraphErrors.h>
64 #include "TTreeStream.h"
65 #include "TGeoManager.h"
67 #include "AliAnalysisManager.h"
68 #include "AliTrackReference.h"
69 #include "AliTrackPointArray.h"
70 #include "AliCDBManager.h"
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"
82 #include "AliTRDtrackInfo/AliTRDclusterInfo.h"
83 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
84 #include "AliTRDtrackingResolution.h"
86 ClassImp(AliTRDtrackingResolution)
88 //________________________________________________________
89 AliTRDtrackingResolution::AliTRDtrackingResolution()
90 :AliTRDrecoTask("Resolution", "Tracking Resolution")
98 ,fTrkltPhiResiduals(0x0)
100 ,fTrkltResolution(0x0)
102 fReconstructor = new AliTRDReconstructor();
103 fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
104 fGeo = new AliTRDgeometry();
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());
115 //________________________________________________________
116 AliTRDtrackingResolution::~AliTRDtrackingResolution()
118 if(fGraphS){fGraphS->Delete(); delete fGraphS;}
119 if(fGraphM){fGraphM->Delete(); delete fGraphM;}
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;}
127 fClResolution->Delete();
128 delete fClResolution;
130 if(fTrkltResolution){fTrkltResolution->Delete(); delete fTrkltResolution;}
134 //________________________________________________________
135 void AliTRDtrackingResolution::CreateOutputObjects()
137 // spatial resolution
138 OpenFile(0, "RECREATE");
140 fContainer = Histos();
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);
154 //________________________________________________________
155 void AliTRDtrackingResolution::Exec(Option_t *opt)
157 fClResiduals->Delete();
158 fTrkltResiduals->Delete();
159 fTrkltPhiResiduals->Delete();
160 fClResolution->Delete();
161 fTrkltResolution->Delete();
163 AliTRDrecoTask::Exec(opt);
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);
172 //________________________________________________________
173 TH1* AliTRDtrackingResolution::PlotClusterResiduals(const AliTRDtrackV1 *track)
175 if(track) fTrack = track;
177 AliWarning("No Track defined.");
181 if(!(h = ((TH2I*)fContainer->At(kClusterResidual)))){
182 AliWarning("No output histogram defined.");
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();
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));
210 //dy = trklt.GetYat(c->GetX()) - c->GetY();
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;
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
240 //________________________________________________________
241 TH1* AliTRDtrackingResolution::PlotTrackletResiduals(const AliTRDtrackV1 *track)
243 if(track) fTrack = track;
245 AliWarning("No Track defined.");
249 if(!(h = ((TH2I*)fContainer->At(kTrackletYResidual)))){
250 AliWarning("No output histogram defined.");
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));
263 // AliRieman fRim(fTrack->GetNumberOfClusters());
264 // Float_t x[AliTRDgeometry::kNlayer] = {-1., -1., -1., -1., -1., -1.}, y[AliTRDgeometry::kNlayer], dydx[AliTRDgeometry::kNlayer];
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);
276 // tracklet->Fit(kTRUE);
278 // x[il] = tracklet->GetX0();
279 // y[il] = tracklet->GetYfit(0)-tracklet->GetYfit(1)*(tracklet->GetX0()-x[il]);
280 // dydx[il] = tracklet->GetYref(1);
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);
290 // if(fDebugLevel>=1){
291 // Double_t sigmay = fRim.GetErrY(x[il]);
292 // Float_t yt = fRim.GetYat(x[il]);
293 // (*fDebugStream) << "TrkltResiduals"
298 // << "dydx=" << dydx[il]
300 // << "sigmay=" << sigmay
306 //________________________________________________________
307 TH1* AliTRDtrackingResolution::PlotTrackletPhiResiduals(const AliTRDtrackV1 *track)
309 if(track) fTrack = track;
311 AliWarning("No Track defined.");
315 if(!(h = ((TH2I*)fContainer->At(kTrackletPhiResidual)))){
316 AliWarning("No output histogram defined.");
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)));
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);
346 // tracklet->Fit(kTRUE);
348 // x[il] = tracklet->GetX0();
349 // y[il] = tracklet->GetYfit(0)-tracklet->GetYfit(1)*(tracklet->GetX0()-x[il]);
350 // dydx[il] = tracklet->GetYref(1);
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;
359 // dydx_ref = fRim.GetDYat(x[il]);
360 // dydx_fit = tracklet->GetYfit(1);
362 // phiref = TMath::ATan(dydx_ref);//*TMath::RadToDeg();
363 // phifit = TMath::ATan(dydx_fit);//*TMath::RadToDeg();
365 // phidiff = phiref-phifit; /*/sigma_phi*/;
367 // h->Fill(dydx_ref, phidiff);
370 // if(fDebugLevel>=1){
371 // (*fDebugStream) << "TrkltPhiResiduals"
373 // << "dydx_fit=" << dydx_fit
374 // << "dydx_ref=" << dydx_ref
375 // << "phiref=" << phiref
376 // << "phifit=" << phifit
377 // << "phidiff=" << phidiff
384 //________________________________________________________
385 TH1* AliTRDtrackingResolution::PlotResolution(const AliTRDtrackV1 *track)
388 AliWarning("No MC defined. Results will not be available.");
391 if(track) fTrack = track;
393 AliWarning("No Track defined.");
397 if(!(h = ((TH2I*)fContainer->At(kClusterResolution)))){
398 AliWarning("No output histogram defined.");
401 if(!(h = ((TH2I*)fContainer->At(kTrackletYResolution)))){
402 AliWarning("No output histogram defined.");
405 if(!(h = ((TH2I*)fContainer->At(kTrackletZResolution)))){
406 AliWarning("No output histogram defined.");
409 if(!(h = ((TH2I*)fContainer->At(kTrackletAngleResolution)))){
410 AliWarning("No output histogram defined.");
413 //printf("called for %d tracklets ... \n", fTrack->GetNumberOfTracklets());
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;
428 // recalculate tracklet based on the MC info
429 AliTRDseedV1 tt(*fTracklet);
432 if(!tt.Fit(kTRUE)) continue;
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);
439 Float_t dphi = TMath::ATan(tt.GetYfit(1)) - TMath::ATan(dydx);
441 Bool_t cross = fTracklet->GetNChange();
443 Double_t *xyz = tt.GetCrossXYZ();
444 dz = xyz[2] - (z0 - (x0 - xyz[0])*dzdx) ;
445 ((TH2I*)fContainer->At(kTrackletZResolution))->Fill(dzdx, dz);
449 ((TH2I*)fContainer->At(kTrackletYResolution))->Fill(dydx, dy);
450 ((TH2I*)fContainer->At(kTrackletAngleResolution))->Fill(dydx, dphi*TMath::RadToDeg());
454 Float_t p = fMC->GetTrackRef() ? fMC->GetTrackRef()->P() : -1.;
455 (*fDebugStream) << "TrkltResolution"
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();
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();
483 Float_t yt = y0 - dx*dydx;
484 Float_t zt = z0 - dx*dzdx;
485 dy = yt - (yc - tilt*(zc-zt));
488 if(q>20. && q<250.) ((TH2I*)fContainer->At(kClusterResolution))->Fill(dydx, dy);
492 Float_t d = zr0 - zt;
493 d -= ((Int_t)(2 * d)) / 2.0;
494 if (d > 0.25) d = 0.5 - d;
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);
506 (*fDebugStream) << "ClusterResolution"
507 <<"clInfo.=" << clInfo
516 //________________________________________________________
517 void AliTRDtrackingResolution::GetRefFigure(Int_t ifig)
521 TGraphErrors *g = 0x0;
523 case kClusterResidual:
524 if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
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;
533 b = new TBox(-.15, -.5, .15, 1.);
534 b->SetFillStyle(3002);b->SetFillColor(kGreen);
535 b->SetLineColor(0); b->Draw();
537 case kTrackletYResidual:
538 if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
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;
547 b = new TBox(-.15, -.5, .15, 1.);
548 b->SetFillStyle(3002);b->SetFillColor(kGreen);
549 b->SetLineColor(0); b->Draw();
551 case kTrackletPhiResidual:
552 if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
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;
561 b = new TBox(-.15, -.5, .15, 1.);
562 b->SetFillStyle(3002);b->SetFillColor(kGreen);
563 b->SetLineColor(0); b->Draw();
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)");
573 if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
575 b = new TBox(-.15, -.5, .15, 1.);
576 b->SetFillStyle(3002);b->SetFillColor(kGreen);
577 b->SetLineColor(0); b->Draw();
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)");
587 if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
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)");
598 if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
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)");
609 if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
613 AliInfo(Form("Reference plot [%d] not implemented yet", ifig));
616 AliInfo(Form("Reference plot [%d] missing result", ifig));
620 //________________________________________________________
621 Bool_t AliTRDtrackingResolution::PostProcess()
623 //fContainer = dynamic_cast<TObjArray*>(GetOutputData(0));
625 Printf("ERROR: list not available");
628 fNRefFigures = fContainer->GetEntriesFast();
630 fGraphS = new TObjArray(fNRefFigures);
634 fGraphM = new TObjArray(fNRefFigures);
640 TGraphErrors *gm = 0x0, *gs = 0x0;
643 TF1 f("f1", "gaus", -.5, .5);
645 TF1 fb("fb", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]", -.5, .5);
647 TF1 fc("fc", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]*exp(-0.5*((x-[4])/[5])**2)", -.5, .5);
650 if(IsVisual()) c = new TCanvas("c", Form("%s Visual", GetName()), 500, 500);
652 sprintf(opt, "%c%c", IsVerbose() ? ' ' : 'Q', IsVisual() ? ' ': 'N');
655 //PROCESS RESIDUAL DISTRIBUTIONS
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;
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);}
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));
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;
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);}
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));
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;
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);}
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));
751 //PROCESS RESOLUTION DISTRIBUTIONS
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;
773 if(IsVisual()){c->cd(); c->SetLogy();}
774 h->Fit(&f, opt, "", -0.5, 0.5);
776 printf("phi[%d] mean[%e] sigma[%e]\n\n", iphi, 10.*f.GetParameter(1), 10.*f.GetParameter(2));
778 if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
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));
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;
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);}
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));
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;
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);}
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));
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;
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);}
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));
886 //________________________________________________________
887 void AliTRDtrackingResolution::Terminate(Option_t *)
894 if(HasPostProcess()) PostProcess();
897 //________________________________________________________
898 void AliTRDtrackingResolution::AdjustF1(TH1 *h, TF1 *f)
900 // Helper function to avoid duplication of code
901 // Make first guesses on the fit parameters
903 // find the intial parameters of the fit !! (thanks George)
904 Int_t nbinsy = Int_t(.5*h->GetNbinsX());
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);
910 f->SetParLimits(1, -.2, .2);
911 f->SetParameter(1, -0.1);
913 f->SetParLimits(2, 0., 4.e-1);
914 f->SetParameter(2, 2.e-2);
915 if(f->GetNpar() <= 4) return;
917 f->SetParLimits(3, 0., sum);
918 f->SetParameter(3, .1*sum);
920 f->SetParLimits(4, -.3, .3);
921 f->SetParameter(4, 0.);
923 f->SetParLimits(5, 0., 1.e2);
924 f->SetParameter(5, 2.e-1);
927 //________________________________________________________
928 TObjArray* AliTRDtrackingResolution::Histos()
930 if(fContainer) return fContainer;
932 fContainer = new TObjArray(7);
933 //fContainer->SetOwner(kTRUE);
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");
943 fContainer->AddAt(h, kClusterResidual);
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");
952 fContainer->AddAt(h, kTrackletYResidual);
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");
961 fContainer->AddAt(h, kTrackletPhiResidual);
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");
973 fContainer->AddAt(h, kClusterResolution);
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");
982 fContainer->AddAt(h, kTrackletYResolution);
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");
991 fContainer->AddAt(h, kTrackletZResolution);
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");
1000 fContainer->AddAt(h, kTrackletAngleResolution);
1006 //________________________________________________________
1007 void AliTRDtrackingResolution::SetRecoParam(AliTRDrecoParam *r)
1010 fReconstructor->SetRecoParam(r);