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 // Reconstruction QA //
23 // Markus Fasel <M.Fasel@gsi.de> //
25 ////////////////////////////////////////////////////////////////////////////
30 #include <TObjArray.h>
35 #include <TGraphErrors.h>
37 #include "TTreeStream.h"
38 #include "TGeoManager.h"
40 #include "AliAnalysisManager.h"
41 #include "AliTrackReference.h"
42 #include "AliTrackPointArray.h"
43 #include "AliCDBManager.h"
45 #include "AliTRDcluster.h"
46 #include "AliTRDseedV1.h"
47 #include "AliTRDtrackV1.h"
48 #include "AliTRDtrackerV1.h"
49 #include "AliTRDReconstructor.h"
50 #include "AliTRDrecoParam.h"
52 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
53 #include "AliTRDtrackingResolution.h"
55 ClassImp(AliTRDtrackingResolution)
57 //________________________________________________________
58 AliTRDtrackingResolution::AliTRDtrackingResolution()
59 :AliTRDrecoTask("Resolution", "Tracking Resolution")
62 fReconstructor = new AliTRDReconstructor();
63 fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
66 //________________________________________________________
67 AliTRDtrackingResolution::~AliTRDtrackingResolution()
69 delete fReconstructor;
70 if(gGeoManager) delete gGeoManager;
74 //________________________________________________________
75 void AliTRDtrackingResolution::CreateOutputObjects()
78 OpenFile(0, "RECREATE");
80 fContainer = Histos();
82 // cluster to tracklet residuals [2]
83 fContainer->AddAt(new TH2I("fYClRes", "Clusters Residuals", 21, -21., 21., 100, -.5, .5), kClusterYResidual);
84 // tracklet to Riemann fit residuals [2]
85 fContainer->AddAt(new TH2I("fYTrkltRRes", "Tracklet Riemann Residuals", 21, -21., 21., 100, -.5, .5), kTrackletRiemanYResidual);
86 fContainer->AddAt(new TH2I("fAngleTrkltRRes", "Tracklet Riemann Angular Residuals", 21, -21., 21., 100, -.5, .5), kTrackletRiemanAngleResidual);
87 fContainer->AddAt(new TH2I("fYTrkltKRes", "Tracklet Kalman Residuals", 21, -21., 21., 100, -.5, .5), kTrackletKalmanYResidual);
88 fContainer->AddAt(new TH2I("fAngleTrkltKRes", "Tracklet Kalman Angular Residuals", 21, -21., 21., 100, -.5, .5), kTrackletKalmanAngleResidual);
92 // tracklet resolution [0]
93 fContainer->AddAt(new TH2I("fY", "Tracklet Resolution", 21, -21., 21., 100, -.5, .5), kTrackletYResolution);
94 // tracklet angular resolution [1]
95 fContainer->AddAt(new TH2I("fPhi", "Tracklet Angular Resolution", 21, -21., 21., 100, -10., 10.), kTrackletAngleResolution);
97 // Riemann track resolution [y, z, angular]
98 fContainer->AddAt(new TH2I("fYRT", "Track Riemann Y Resolution", 21, -21., 21., 100, -.5, .5), kTrackRYResolution);
99 fContainer->AddAt(new TH2I("fZRT", "Track Riemann Z Resolution", 21, -21., 21., 100, -.5, .5), kTrackRZResolution);
100 fContainer->AddAt(new TH2I("fPhiRT", "Track Riemann Angular Resolution", 21, -21., 21., 100, -10., 10.), kTrackRAngleResolution);
102 // Kalman track resolution [y, z, angular]
103 fContainer->AddAt(new TH2I("fYKT", "", 21, -21., 21., 100, -.5, .5), kTrackKYResolution);
104 fContainer->AddAt(new TH2I("fZKT", "", 21, -21., 21., 100, -.5, .5), kTrackKZResolution);
105 fContainer->AddAt(new TH2I("fPhiKT", "", 21, -21., 21., 100, -10., 10.), kTrackKAngleResolution);
108 // CREATE GRAPHS for DISPLAY
110 // define iterator over graphs
111 Int_t jgraph = (Int_t)kGraphStart;
112 TH2I *h2 = (TH2I *)(fContainer->At(kClusterYResidual));
113 // clusters tracklet residuals (mean-phi)
114 TH1 *h = new TH1I("h", "", 100, -40., 40.);
115 h->GetXaxis()->SetTitle("#Phi [deg]");
116 h->GetYaxis()->SetTitle("Clusters Residuals : #sigma/#mu [mm]");
117 h->GetYaxis()->SetRangeUser(-.05, 1.);
118 fContainer->AddAt(h, jgraph++);
120 TGraphErrors *g = new TGraphErrors(h2->GetNbinsX());
121 g->SetLineColor(kGreen);
122 g->SetMarkerStyle(22);
123 g->SetMarkerColor(kGreen);
124 g->SetNameTitle("clm", "Residuals Clusters-Tracklet Mean");
125 fContainer->AddAt(g, jgraph++);
127 // clusters tracklet residuals (sigma-phi)
128 g = new TGraphErrors(h2->GetNbinsX());
129 g->SetLineColor(kRed);
130 g->SetMarkerStyle(23);
131 g->SetMarkerColor(kRed);
132 g->SetNameTitle("cls", "Residuals Clusters-Tracklet Sigma");
133 fContainer->AddAt(g, jgraph++);
136 // tracklet y resolution
137 h2 = (TH2I*)fContainer->At(kTrackletYResolution);
138 h = new TH1I("h", "", 100, -40., 40.);
139 h->GetXaxis()->SetTitle("#Phi [deg]");
140 h->GetYaxis()->SetTitle("Tracklet Resolution : #sigma/#mu [mm]");
141 h->GetYaxis()->SetRangeUser(-.05, 1.);
142 fContainer->AddAt(h, jgraph++);
144 g = new TGraphErrors(h2->GetNbinsX());
145 g->SetLineColor(kGreen);
146 g->SetMarkerStyle(22);
147 g->SetMarkerColor(kGreen);
148 g->SetNameTitle("trkltym", "Resolution Tracklet Y Mean");
149 fContainer->AddAt(g, jgraph++);
150 g = new TGraphErrors(h2->GetNbinsX());
151 g->SetLineColor(kRed);
152 g->SetMarkerStyle(22);
153 g->SetMarkerColor(kRed);
154 g->SetNameTitle("trkltys", "Resolution Tracklet Y Sigma");
155 fContainer->AddAt(g, jgraph++);
157 // tracklet phi resolution
158 h2 = (TH2I*)fContainer->At(kTrackletAngleResolution);
159 h = new TH1I("h", "", 100, -40., 40.);
160 h->GetXaxis()->SetTitle("#Phi [deg]");
161 h->GetYaxis()->SetTitle("Tracklet Angular Resolution : #sigma/#mu [deg]");
162 h->GetYaxis()->SetRangeUser(-.05, .2);
163 fContainer->AddAt(h, jgraph++);
165 g = new TGraphErrors(h2->GetNbinsX());
166 g->SetLineColor(kGreen);
167 g->SetMarkerStyle(22);
168 g->SetMarkerColor(kGreen);
169 g->SetNameTitle("trkltam", "Resolution Tracklet Y Mean");
170 fContainer->AddAt(g, jgraph++);
171 g = new TGraphErrors(h2->GetNbinsX());
172 g->SetLineColor(kRed);
173 g->SetMarkerStyle(22);
174 g->SetMarkerColor(kRed);
175 g->SetNameTitle("trkltas", "Angle Resolution Tracklet Sigma");
176 fContainer->AddAt(g, jgraph++);
180 //________________________________________________________
181 void AliTRDtrackingResolution::Exec(Option_t *)
183 // spatial Resolution: res = pos_{Tracklet}(x = x_{Anode wire}) - pos_{TrackRef}(x = x_{Anode wire})
184 // angular Resolution: res = Tracklet angle - TrackRef Angle
186 Int_t nTrackInfos = fTracks->GetEntriesFast();
187 if(fDebugLevel>=2) printf("Number of Histograms: %d\n", Histos()->GetEntries());
189 Double_t p, dy, dphi, dymc, dzmc, dphimc;
190 Float_t fP[kNLayers], fY[kNLayers], fZ[kNLayers], fPhi[kNLayers], fTheta[kNLayers]; // phi/theta angle per layer
191 Bool_t fMCMap[kNLayers], fLayerMap[kNLayers]; // layer map
193 AliTrackPoint tr[kNLayers], tk[kNLayers];
194 AliExternalTrackParam *fOp = 0x0;
195 AliTRDtrackV1 *fTrack = 0x0;
196 AliTRDtrackInfo *fInfo = 0x0;
197 if(fDebugLevel>=2) printf("Number of TrackInfos: %d\n", nTrackInfos);
198 for(Int_t iTI = 0; iTI < nTrackInfos; iTI++){
199 // check if ESD and MC-Information are available
200 if(!(fInfo = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(iTI)))) continue;
201 if(!(fTrack = fInfo->GetTRDtrack())) continue;
202 if(!(fOp = fInfo->GetOuterParam())) continue;
204 if(fDebugLevel>=3) printf("\tDoing track[%d] NTrackRefs[%d]\n", iTI, fInfo->GetNTrackRefs());
208 memset(fP, 0, kNLayers*sizeof(Float_t));
209 memset(fY, 0, kNLayers*sizeof(Float_t));
210 memset(fZ, 0, kNLayers*sizeof(Float_t));
211 memset(fPhi, 0, kNLayers*sizeof(Float_t));
212 memset(fTheta, 0, kNLayers*sizeof(Float_t));
213 memset(fLayerMap, 0, kNLayers*sizeof(Bool_t));
214 memset(fMCMap, 0, kNLayers*sizeof(Bool_t));
215 AliTRDseedV1 *fTracklet = 0x0;
216 for(Int_t iplane = 0; iplane < kNLayers; iplane++){
217 if(!(fTracklet = fTrack->GetTracklet(iplane))) continue;
218 if(!fTracklet->IsOK()) continue;
221 fLayerMap[iplane] = kTRUE;
222 tr[npts].SetXYZ(fTracklet->GetX0(), 0., 0.);
223 tk[npts].SetXYZ(fTracklet->GetX0(), fTracklet->GetYfit(0), fTracklet->GetZfit(0));
226 if(fDebugLevel>=4) printf("\t\tLy[%d] X0[%6.3f] Ncl[%d]\n", iplane, fTracklet->GetX0(), fTracklet->GetN());
228 // define reference values
230 fY[iplane] = fTracklet->GetYref(0);
231 fZ[iplane] = fTracklet->GetZref(0);
232 fPhi[iplane] = TMath::ATan(fTracklet->GetYref(1));
233 fTheta[iplane] = TMath::ATan(fTracklet->GetZref(1));
236 // RESOLUTION (compare to MC)
238 if(fInfo->GetNTrackRefs() >= 2){
239 Double_t pmc, ymc, zmc, phiMC, thetaMC;
240 if(Resolution(fTracklet, fInfo, pmc, ymc, zmc, phiMC, thetaMC)){
241 fMCMap[iplane] = kTRUE;
245 fPhi[iplane] = phiMC;
246 fTheta[iplane] = thetaMC;
250 Float_t phi = fPhi[iplane]*TMath::RadToDeg();
251 Float_t theta = fTheta[iplane]*TMath::RadToDeg();
253 // Do clusters residuals
254 if(!fTracklet->Fit(kFALSE)) continue;
255 AliTRDcluster *c = 0x0;
256 for(Int_t ic=AliTRDseed::knTimebins-1; ic>=0; ic--){
257 if(!(c = fTracklet->GetClusters(ic))) continue;
259 dy = fTracklet->GetYat(c->GetX()) - c->GetY();
260 ((TH2I*)fContainer->At(kClusterYResidual))->Fill(phi, dy);
262 Float_t q = c->GetQ();
263 (*fDebugStream) << "ResidualClusters"
275 // this protection we might drop TODO
276 if(fTrack->GetNumberOfTracklets() < 6) continue;
278 AliTRDtrackerV1::FitRiemanTilt(fTrack, 0x0, kTRUE, npts, tr);
280 for(Int_t ip=0; ip<kNLayers; ip++){
281 if(!fLayerMap[ip]) continue;
282 fTracklet = fTrack->GetTracklet(ip);
283 // recalculate fit based on the new tilt correction
286 dy = fTracklet->GetYfit(0) - tr[iref].GetY();
287 ((TH2I*)fContainer->At(kTrackletRiemanYResidual))->Fill(fPhi[ip]*TMath::RadToDeg(), dy);
289 dphi = fTracklet->GetYfit(1)- fTracklet->GetYref(1);
290 ((TH2I*)fContainer->At(kTrackletRiemanAngleResidual))->Fill(fPhi[ip]*TMath::RadToDeg(), dphi);
293 dymc = fY[ip] - tr[iref].GetY();
294 ((TH2I*)fContainer->At(kTrackRYResolution))->Fill(fPhi[ip]*TMath::RadToDeg(), dymc);
296 dzmc = fZ[ip] - tr[iref].GetZ();
297 ((TH2I*)fContainer->At(kTrackRZResolution))->Fill(fPhi[ip]*TMath::RadToDeg(), dzmc);
299 dphimc = fPhi[ip] - fTracklet->GetYfit(1);
300 ((TH2I*)fContainer->At(kTrackRAngleResolution))->Fill(fPhi[ip]*TMath::RadToDeg(), dphimc);
306 (*fDebugStream) << "RiemannTrack"
308 << "mc=" << fMCMap[ip]
310 << "phi=" << fPhi[ip]
311 << "tht=" << fTheta[ip]
316 << "dphimc="<< dphimc
321 // if(!gGeoManager) TGeoManager::Import("geometry.root");
322 // AliTRDtrackerV1::FitKalman(fTrack, 0x0, kFALSE, nc, tr);
323 // for(Int_t ip=0; ip<nc; ip++){
324 // dy = cl[ip].GetY() - tr[ip].GetY();
325 // ((TH2I*)fContainer->At(kTrackletKalmanYResidual))->Fill(phi*TMath::RadToDeg(), dy);
326 // dz = cl[ip].GetZ() - tr[ip].GetZ();
327 // if(fDebugLevel>=1){
328 // (*fDebugStream) << "KalmanTrack"
331 // /* << "phi=" << phi
332 // << "theta=" << theta
333 // << "dphi=" << dphi*/
340 PostData(0, fContainer);
343 //________________________________________________________
344 void AliTRDtrackingResolution::GetRefFigure(Int_t ifig, Int_t &first, Int_t &last, Option_t */*opt*/)
348 first = (Int_t)kGraphStart; last = first+3;
351 first = (Int_t)kGraphStart+3; last = first+3;
354 first = (Int_t)kGraphStart+6; last = first+3;
357 first = (Int_t)kGraphStart; last = first;
363 //________________________________________________________
364 Bool_t AliTRDtrackingResolution::Resolution(AliTRDseedV1 *tracklet, AliTRDtrackInfo *fInfo, Double_t &p, Double_t &ymc, Double_t &zmc, Double_t &phi, Double_t &theta)
367 AliTrackReference *fTrackRefs[2] = {0x0, 0x0}, *tempTrackRef = 0x0;
369 // check for 2 track ref where the radial position has a distance less than 3.7mm
371 for(Int_t itr = 0; itr < fInfo->GetNTrackRefs(); itr++){
372 if(!(tempTrackRef = fInfo->GetTrackRef(itr))) continue;
373 if(fDebugLevel>=5) printf("TrackRef %d: x = %f\n", itr, tempTrackRef->LocalX());
374 if(TMath::Abs(tracklet->GetX0() - tempTrackRef->LocalX()) > 3.7) continue;
375 fTrackRefs[nFound++] = tempTrackRef;
376 if(nFound == 2) break;
379 if(fDebugLevel>=4) printf("\t\tFound track crossing [%d] refX[%6.3f]\n", nFound, nFound>0 ? fTrackRefs[0]->LocalX() : 0.);
382 // We found 2 track refs for the tracklet, get y and z at the anode wire by a linear approximation
386 Double_t dx = fTrackRefs[1]->LocalX() - fTrackRefs[0]->LocalX();
388 if(fDebugLevel>=4) printf("\t\ttrack ref in the wrong order refX0[%6.3f] refX1[%6.3f]\n", fTrackRefs[0]->LocalX(), fTrackRefs[1]->LocalX());
391 Double_t dydx = (fTrackRefs[1]->LocalY() - fTrackRefs[0]->LocalY()) / dx;
392 Double_t dzdx = (fTrackRefs[1]->Z() - fTrackRefs[0]->Z()) / dx;
393 Double_t dx0 = fTrackRefs[1]->LocalX() - tracklet->GetX0();
394 ymc = fTrackRefs[1]->LocalY() - dydx*dx0;
395 zmc = fTrackRefs[1]->Z() - dzdx*dx0;
397 // recalculate tracklet based on the MC info
398 tracklet->SetZref(0, zmc);
399 tracklet->SetZref(1, -dzdx); // TODO
401 Double_t dy = tracklet->GetYfit(0) - ymc;
402 Double_t dz = tracklet->GetZfit(0) - zmc;
404 //res_y *= 100; // in mm
405 p = fTrackRefs[0]->P();
407 phi = TMath::ATan(dydx);
408 theta = TMath::ATan(dzdx);
409 Double_t dphi = TMath::ATan(tracklet->GetYfit(1)) - phi;
410 if(fDebugLevel>=4) printf("\t\tdx[%6.4f] dy[%6.4f] dz[%6.4f] dphi[%6.4f] \n", dx, dy, dz, dphi);
413 if(TMath::Abs(dx-3.7)<1.E-3){
414 ((TH2I*)fContainer->At(kTrackletYResolution))->Fill(phi*TMath::RadToDeg(), dy);
415 ((TH2I*)fContainer->At(kTrackletAngleResolution))->Fill(phi*TMath::RadToDeg(), dphi*TMath::RadToDeg());
419 Int_t iplane = tracklet->GetPlane();
420 (*fDebugStream) << "ResolutionTrklt"
431 << "tracklet.="<<tracklet
438 //________________________________________________________
439 Bool_t AliTRDtrackingResolution::PostProcess()
441 //fContainer = dynamic_cast<TObjArray*>(GetOutputData(0));
444 Printf("ERROR: list not available");
450 TGraphErrors *gm = 0x0, *gs = 0x0;
451 TF1 f("f1", "gaus", -.5, .5);
452 // define iterator over graphs
453 Int_t jgraph = (Int_t)kGraphStart;
455 //PROCESS RESIDUAL DISTRIBUTIONS
457 // Clusters residuals
458 h2 = (TH2I *)(fContainer->At(kClusterYResidual));
459 jgraph++; //skip the frame histo
460 gm = (TGraphErrors*)fContainer->At(jgraph++);
461 gs = (TGraphErrors*)fContainer->At(jgraph++);
462 for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
463 Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
464 Double_t dphi = h2->GetXaxis()->GetBinWidth(ibin)/2;
465 h = h2->ProjectionY("py", ibin, ibin);
466 h->Fit(&f, "QN", "", -0.5, 0.5);
467 gm->SetPoint(ibin - 1, phi, 10.*f.GetParameter(1));
468 gm->SetPointError(ibin - 1, dphi, 10.*f.GetParError(1));
469 gs->SetPoint(ibin - 1, phi, 10.*f.GetParameter(2));
470 gs->SetPointError(ibin - 1, dphi, 10.*f.GetParError(2));
475 //PROCESS RESOLUTION DISTRIBUTIONS
477 // tracklet y resolution
478 h2 = (TH2I*)fContainer->At(kTrackletYResolution);
479 jgraph++; //skip the frame histo
480 gm = (TGraphErrors*)fContainer->At(jgraph++);
481 gs = (TGraphErrors*)fContainer->At(jgraph++);
482 for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
483 Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
484 f.SetParameter(1, 0.);f.SetParameter(2, 2.e-2);
485 h = h2->ProjectionY("py", iphi, iphi);
486 h->Fit(&f, "QN", "", -.5, .5);
487 Int_t jphi = iphi -1;
488 gm->SetPoint(jphi, phi, 10.*f.GetParameter(1));
489 gm->SetPointError(jphi, 0., 10.*f.GetParError(1));
490 gs->SetPoint(jphi, phi, 10.*f.GetParameter(2));
491 gs->SetPointError(jphi, 0., 10.*f.GetParError(2));
495 // tracklet phi resolution
496 h2 = (TH2I*)fContainer->At(kTrackletAngleResolution);
497 jgraph++; //skip the frame histo
498 gm = (TGraphErrors*)fContainer->At(jgraph++);
499 gs = (TGraphErrors*)fContainer->At(jgraph++);
500 for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
501 Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
502 f.SetParameter(1, 0.);f.SetParameter(2, 2.e-2);
503 h = h2->ProjectionY("py", iphi, iphi);
504 h->Fit(&f, "QN", "", -.5, .5);
505 Int_t jphi = iphi -1;
506 gm->SetPoint(jphi, phi, f.GetParameter(1));
507 gm->SetPointError(jphi, 0., f.GetParError(1));
508 gs->SetPoint(jphi, phi, f.GetParameter(2));
509 gs->SetPointError(jphi, 0., f.GetParError(2));
518 //________________________________________________________
519 void AliTRDtrackingResolution::Terminate(Option_t *)
526 if(HasPostProcess()) PostProcess();
529 //________________________________________________________
530 TObjArray* AliTRDtrackingResolution::Histos()
532 if(!fContainer) fContainer = new TObjArray(25);
537 //________________________________________________________
538 void AliTRDtrackingResolution::SetRecoParam(AliTRDrecoParam *r)
541 fReconstructor->SetRecoParam(r);