]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG1/TRD/AliTRDclusterResolution.cxx
initialize magnetic field in infoGen (moved from clusterResolution)
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDclusterResolution.cxx
CommitLineData
1ee39b3a 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: AliTRDclusterResolution.cxx */
17
18///////////////////////////////////////////////////////////////////////////////
19// //
20// TRD cluster error parameterization //
21// //
22// This class is designed to produce the reference plots for a detailed study//
23// and parameterization of TRD cluster errors. The following effects are taken//
24// into account : //
25// - dependence with the total charge of the cluster //
26// - dependence with the distance from the center pad. This is monitored
27// for each layer individually since the pad size varies with layer
28// - dependence with the drift length - here the influence of anisochronity
29// and diffusion are searched
30// - dependence with the distance to the anode wire - anisochronity effects
31// - dependence with track angle (for y resolution)
32// The correlation between effects is taken into account.
33//
34// Since magnetic field plays a very important role in the TRD measurement
35// the ExB correction is forced by the setter function SetExB(Int_t). The
36// argument is the detector index, if none is specified all will be
37// considered.
38//
39// Two cases are of big importance.
40// - comparison with MC
41// - comparison with Kalman fit. In this case the covariance matrix of the
42// Kalman fit are needed.
43//
44// The functionalities implemented in this class are based on the storage
45// class AliTRDclusterInfo.
46//
47// The Method
48// ----------
49//
50// The method to disentangle s_y and s_x is based on the relation (see also fig.)
51// BEGIN_LATEX
52// #sigma^{2} = #sigma^{2}_{y} + tg^{2}(#alpha_{L})*#sigma^{2}_{x_{d}} + tg^{2}(#phi-#alpha_{L})*(#sigma^{2}_{x_{d}}+#sigma^{2}_{x_{c}})
53// END_LATEX
54// with
55// BEGIN_LATEX
56// #sigma^{2}_{x_{c}} #approx 0
57// END_LATEX
58// we suppose the chamber is well calibrated for t_{0} and aligned in
59// radial direction.
60//
61// Clusters can be radially shifted due to three causes:
62// - globally shifted - due to residual misalignment/miscalibration(t0)
63// - locally shifted - due to different local drift velocity from the mean
64// - randomly shifted - due to neighboring (radial direction) clusters
65// charge induced by asymmetry of the TRF.
66//
67// We estimate this effects by the relations:
68// BEGIN_LATEX
69// #mu_{y} = tg(#alpha_{L})*#Delta x_{d}(...) + tg(#phi-#alpha_{L})*(#Delta x_{c}(...) + #Delta x_{d}(...))
70// END_LATEX
71// where
72// BEGIN_LATEX
73// #Delta x_{d}(...) = (<v_{d}> + #delta v_{d}(x_{d}, d)) * (t + t^{*}(Q))
74// END_LATEX
75// and we specified explicitely the variation of drift velocity parallel
76// with the track (x_{d}) and perpendicular to it due to anisochronity (d).
77//
78// For estimating the contribution from asymmetry of TRF the following
79// parameterization is being used
80// BEGIN_LATEX
81// t^{*}(Q) = #delta_{0} * #frac{Q_{t+1} - Q_{t-1}}{Q_{t-1} + Q_{t} + Q_{t+1}}
82// END_LATEX
83//
84//
85// Clusters can also be r-phi shifted due to:
86// - wrong PRF or wrong cuts at digits level
87//The following correction is applied :
88// BEGIN_LATEX
89// <#Delta y> = a + b * sin(c*y_{pw})
90// END_LATEX
91
92// The Models
93//
94// Parameterization against total charge
95//
96// Obtained for B=0T at phi=0. All other effects integrated out.
97// BEGIN_LATEX
98// #sigma^{2}_{y}(Q) = #sigma^{2}_{y}(...) + b(#frac{1}{Q} - #frac{1}{Q_{0}})
99// END_LATEX
100// For B diff 0T the error of the average ExB correction error has to be subtracted !!
101//
102// Parameterization Sx
103//
104// The parameterization of the error in the x direction can be written as
105// BEGIN_LATEX
106// #sigma_{x} = #sigma_{x}^{||} + #sigma_{x}^{#perp}
107// END_LATEX
108//
109// where the parallel component is given mainly by the TRF width while
110// the perpendicular component by the anisochronity. The model employed for
111// the parallel is gaus(0)+expo(3) with the following parameters
112// 1 C 5.49018e-01 1.23854e+00 3.84540e-04 -8.21084e-06
113// 2 M 7.82999e-01 6.22531e-01 2.71272e-04 -6.88485e-05
114// 3 S 2.74451e-01 1.13815e+00 2.90667e-04 1.13493e-05
115// 4 E1 2.53596e-01 1.08646e+00 9.95591e-05 -2.11625e-05
116// 5 E2 -2.40078e-02 4.26520e-01 4.67153e-05 -2.35392e-04
117//
118// and perpendicular to the track is pol2 with the parameters
119//
120// Par_0 = 0.190676 +/- 0.41785
121// Par_1 = -3.9269 +/- 7.49862
122// Par_2 = 14.7851 +/- 27.8012
123//
124// Parameterization Sy
125//
126// The parameterization of the error in the y direction along track uses
127// BEGIN_LATEX
128// #sigma_{y}^{||} = #sigma_{y}^{0} -a*exp(1/(x-b))
129// END_LATEX
130//
131// with following values for the parameters:
132// 1 sy0 2.60967e-01 2.99652e-03 7.82902e-06 -1.89636e-04
133// 2 a -7.68941e+00 1.87883e+00 3.84539e-04 9.38268e-07
134// 3 b -3.41160e-01 7.72850e-02 1.63231e-05 2.51602e-05
135//
136//==========================================================================
137// Example how to retrive reference plots from the task
138// void steerClErrParam(Int_t fig=0)
139// {
140// gSystem->Load("libANALYSIS.so");
141// gSystem->Load("libTRDqaRec.so");
142//
143// // initialize DB manager
144// AliCDBManager *cdb = AliCDBManager::Instance();
145// cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
146// cdb->SetRun(0);
147// // initialize magnetic field.
148// AliMagFCheb *field=new AliMagFCheb("Maps","Maps", 2, 1., 10., AliMagFCheb::k5kG);
149// AliTracker::SetFieldMap(field, kTRUE);
150//
151// AliTRDclusterResolution *res = new AliTRDclusterResolution();
152// res->SetMCdata();
153// res->Load("TRD.TaskClErrParam.root");
154// res->SetExB();
155// res->SetVisual();
156// //res->SetSaveAs();
157// res->SetProcessCharge(kFALSE);
158// res->SetProcessCenterPad(kFALSE);
159// //res->SetProcessMean(kFALSE);
160// res->SetProcessSigma(kFALSE);
161// if(!res->PostProcess()) return;
162// new TCanvas;
163// res->GetRefFigure(fig);
164// }
165//
166// Authors: //
167// Alexandru Bercuci <A.Bercuci@gsi.de> //
168////////////////////////////////////////////////////////////////////////////
169
170#include "AliTRDclusterResolution.h"
171#include "info/AliTRDclusterInfo.h"
172#include "AliTRDgeometry.h"
801d4d50 173#include "AliTRDpadPlane.h"
1ee39b3a 174#include "AliTRDcluster.h"
5935a6da 175#include "AliTRDseedV1.h"
1ee39b3a 176#include "AliTRDcalibDB.h"
177#include "AliTRDCommonParam.h"
178#include "Cal/AliTRDCalROC.h"
179#include "Cal/AliTRDCalDet.h"
180
801d4d50 181#include "AliESDEvent.h"
1ee39b3a 182#include "AliCDBManager.h"
183
184#include "TROOT.h"
185#include "TObjArray.h"
186#include "TAxis.h"
187#include "TF1.h"
188#include "TLegend.h"
189#include "TGraphErrors.h"
190#include "TLine.h"
191#include "TH2I.h"
192#include "TH3S.h"
193#include "TTree.h"
194#include "TMath.h"
195#include "TLinearFitter.h"
1a68da85 196#include "TGeoGlobalMagField.h"
1ee39b3a 197
198#include "TCanvas.h"
199#include "TSystem.h"
200
201ClassImp(AliTRDclusterResolution)
202
203const Float_t AliTRDclusterResolution::fgkTimeBinLength = 1./ AliTRDCommonParam::Instance()->GetSamplingFrequency();
204//_______________________________________________________
f8f46e4d 205AliTRDclusterResolution::AliTRDclusterResolution()
206 : AliTRDrecoTask()
705f8b0a 207 ,fCanvas(NULL)
208 ,fInfo(NULL)
209 ,fResults(NULL)
56f313bd 210 ,fSubTaskMap(0)
211 ,fUseCalib(7)
f8f46e4d 212 ,fDet(-1)
801d4d50 213 ,fCol(-1)
214 ,fRow(-1)
f8f46e4d 215 ,fExB(0.)
e3147647 216 ,fVdrift(1.5)
5935a6da 217 ,fT0(0.)
e3147647 218 ,fGain(1.)
563d1b38 219 ,fDyRange(1.3)
f8f46e4d 220 ,fLy(0)
5935a6da 221 ,fT(0.)
f8f46e4d 222 ,fX(0.)
223 ,fY(0.)
224 ,fZ(0.)
225{
226// Constructor
705f8b0a 227 SetNameTitle("ClErrCalib", "Cluster Error Parameterization");
563d1b38 228 memset(fR, 0, 4*sizeof(Float_t));
229 memset(fP, 0, 4*sizeof(Float_t));
f8f46e4d 230}
231
705f8b0a 232//_______________________________________________________
233AliTRDclusterResolution::AliTRDclusterResolution(const char *name)
234 : AliTRDrecoTask(name, "Cluster Error Parameterization")
4226db3e 235 ,fCanvas(NULL)
236 ,fInfo(NULL)
237 ,fResults(NULL)
56f313bd 238 ,fSubTaskMap(0)
239 ,fUseCalib(7)
1ee39b3a 240 ,fDet(-1)
801d4d50 241 ,fCol(-1)
242 ,fRow(-1)
1ee39b3a 243 ,fExB(0.)
e3147647 244 ,fVdrift(1.5)
5935a6da 245 ,fT0(0.)
e3147647 246 ,fGain(1.)
563d1b38 247 ,fDyRange(1.3)
1ee39b3a 248 ,fLy(0)
5935a6da 249 ,fT(0.)
1ee39b3a 250 ,fX(0.)
251 ,fY(0.)
252 ,fZ(0.)
253{
254// Constructor
255
256 memset(fR, 0, 4*sizeof(Float_t));
257 memset(fP, 0, 4*sizeof(Float_t));
1ee39b3a 258
259 // By default register all analysis
260 // The user can switch them off in his steering macro
261 SetProcess(kQRes);
262 SetProcess(kCenter);
263 SetProcess(kMean);
264 SetProcess(kSigm);
265}
266
267//_______________________________________________________
268AliTRDclusterResolution::~AliTRDclusterResolution()
269{
270// Destructor
271
272 if(fCanvas) delete fCanvas;
1ee39b3a 273 if(fResults){
274 fResults->Delete();
275 delete fResults;
276 }
277}
278
1ee39b3a 279//_______________________________________________________
f8f46e4d 280void AliTRDclusterResolution::UserCreateOutputObjects()
1ee39b3a 281{
1ee39b3a 282 fContainer = Histos();
068e2c00 283 PostData(1, fContainer);
1ee39b3a 284}
285
286//_______________________________________________________
287Bool_t AliTRDclusterResolution::GetRefFigure(Int_t ifig)
288{
289// Steering function to retrieve performance plots
290
291 if(!fResults) return kFALSE;
4226db3e 292 TLegend *leg = NULL;
293 TList *l = NULL;
294 TObjArray *arr = NULL;
295 TTree *t = NULL;
296 TH2 *h2 = NULL;TH1 *h1 = NULL;
297 TGraphErrors *gm(NULL), *gs(NULL), *gp(NULL);
1ee39b3a 298 switch(ifig){
299 case kQRes:
300 if(!(arr = (TObjArray*)fResults->At(kQRes))) break;
301 if(!(gm = (TGraphErrors*)arr->At(0))) break;
302 if(!(gs = (TGraphErrors*)arr->At(1))) break;
303 if(!(gp = (TGraphErrors*)arr->At(2))) break;
5935a6da 304 leg= new TLegend(.7, .7, .9, .95);
305 leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetFillStyle(0);
306 gs->Draw("apl"); leg->AddEntry(gs, "Sigma / Resolution", "pl");
1ee39b3a 307 gs->GetHistogram()->GetYaxis()->SetRangeUser(-50., 700.);
308 gs->GetHistogram()->SetXTitle("Q [a.u.]");
5935a6da 309 gs->GetHistogram()->SetYTitle("y - x tg(#alpha_{L}) [#mum]");
310 gm->Draw("pl");leg->AddEntry(gm, "Mean / Systematics", "pl");
311 gp->Draw("pl");leg->AddEntry(gp, "Abundance / Probability", "pl");
312 leg->Draw();
1ee39b3a 313 return kTRUE;
314 case kCenter:
315 if(!(arr = (TObjArray*)fResults->At(kCenter))) break;
316 gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
317 ((TVirtualPad*)l->At(0))->cd();
5935a6da 318 ((TTree*)arr->At(0))->Draw(Form("y:t>>h(%d, -0.5, %f, 51, -.51, .51)",AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5),
1ee39b3a 319 "m[0]*(ly==0&&abs(m[0])<1.e-1)", "colz");
320 ((TVirtualPad*)l->At(1))->cd();
321 leg= new TLegend(.7, .7, .9, .95);
322 leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetFillStyle(0);
323 leg->SetHeader("TRD Plane");
324 for(Int_t il = 1; il<=AliTRDgeometry::kNlayer; il++){
325 if(!(gm = (TGraphErrors*)arr->At(il))) return kFALSE;
326 gm->Draw(il>1?"pc":"apc"); leg->AddEntry(gm, Form("%d", il-1), "pl");
327 if(il>1) continue;
5935a6da 328 gm->GetHistogram()->SetXTitle("t_{drift} [tb]");
1ee39b3a 329 gm->GetHistogram()->SetYTitle("#sigma_{y}(x|cen=0) [#mum]");
330 gm->GetHistogram()->GetYaxis()->SetRangeUser(150., 500.);
331 }
332 leg->Draw();
333 return kTRUE;
334 case kSigm:
335 if(!(t = (TTree*)fResults->At(kSigm))) break;
336 t->Draw("z:t>>h2x(23, 0.1, 2.4, 25, 0., 2.5)","sx*(1)", "lego2fb");
337 h2 = (TH2F*)gROOT->FindObject("h2x");
338 printf(" const Double_t sx[24][25]={\n");
339 for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
340 printf(" {");
341 for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
342 printf("%6.4f ", h2->GetBinContent(ix, iy));
343 }
344 printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
345 }
346 printf(" };\n");
347 gPad->Divide(2, 1, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
348 ((TVirtualPad*)l->At(0))->cd();
349 h1 = h2->ProjectionX("hsx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
350 h1->SetYTitle("<#sigma_{x}> [#mum]");
351 h1->SetXTitle("t_{drift} [#mus]");
5935a6da 352 h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); h1->Draw("pc");
1ee39b3a 353
354 t->Draw("z:t>>h2y(23, 0.1, 2.4, 25, 0., 2.5)","sy*(1)", "lego2fb");
355 h2 = (TH2F*)gROOT->FindObject("h2y");
356 printf(" const Double_t sy[24][25]={\n");
357 for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
358 printf(" {");
359 for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
360 printf("%6.4f ", h2->GetBinContent(ix, iy));
361 }
362 printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
363 }
364 printf(" };\n");
365 ((TVirtualPad*)l->At(1))->cd();
366 h1 = h2->ProjectionX("hsy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
367 h1->SetYTitle("<#sigma_{y}> [#mum]");
368 h1->SetXTitle("t_{drift} [#mus]");
5935a6da 369 h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); h1->Draw("pc");
1ee39b3a 370 return kTRUE;
371 case kMean:
372 if(!(t = (TTree*)fResults->At(kMean))) break;
2ba7720d 373 if(!t->Draw(Form("z:t>>h2x(%d, -0.5, %3.1f, %d, 0., 2.5)",
5935a6da 374 AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5, kND),
2ba7720d 375 "dx*(1)", "goff")) break;
1ee39b3a 376 h2 = (TH2F*)gROOT->FindObject("h2x");
5935a6da 377 printf(" const Double_t dx[%d][%d]={\n", AliTRDseedV1::kNtb, kND);
1ee39b3a 378 for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
379 printf(" {");
380 for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
5935a6da 381 printf("%+6.4e, ", h2->GetBinContent(ix, iy));
1ee39b3a 382 }
5935a6da 383 printf("%+6.4e},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
1ee39b3a 384 }
385 printf(" };\n");
5935a6da 386 gPad->Divide(2, 2, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
1ee39b3a 387 ((TVirtualPad*)l->At(0))->cd();
5935a6da 388 h2->Draw("lego2fb");
389 ((TVirtualPad*)l->At(2))->cd();
1ee39b3a 390 h1 = h2->ProjectionX("hdx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
5935a6da 391 h1->SetYTitle("<#deltax> [#mum]");
b9ddd472 392 h1->SetXTitle("t_{drift} [tb]");
5935a6da 393 //h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1);
394 h1->Draw("pc");
1ee39b3a 395
2ba7720d 396 if(!t->Draw(Form("z:t>>h2y(%d, -0.5, %3.1f, %d, 0., 2.5)",
5935a6da 397 AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5, kND),
2ba7720d 398 "dy*(1)", "goff")) break;
1ee39b3a 399 h2 = (TH2F*)gROOT->FindObject("h2y");
5935a6da 400 printf(" const Double_t dy[%d][%d]={\n", AliTRDseedV1::kNtb, kND);
1ee39b3a 401 for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
402 printf(" {");
403 for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
5935a6da 404 printf("%+6.4e ", h2->GetBinContent(ix, iy));
1ee39b3a 405 }
5935a6da 406 printf("%+6.4e},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
1ee39b3a 407 }
408 printf(" };\n");
409 ((TVirtualPad*)l->At(1))->cd();
5935a6da 410 h2->Draw("lego2fb");
411 ((TVirtualPad*)l->At(3))->cd();
1ee39b3a 412 h1 = h2->ProjectionX("hdy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
5935a6da 413 h1->SetYTitle("<#deltay> [#mum]");
b9ddd472 414 h1->SetXTitle("t_{drift} [tb]");
5935a6da 415 //h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1);
416 h1->Draw("pc");
1ee39b3a 417
418 return kTRUE;
419 default:
420 break;
421 }
422 AliWarning("No container/data found.");
423 return kFALSE;
424}
425
426//_______________________________________________________
427TObjArray* AliTRDclusterResolution::Histos()
428{
429// Retrieve histograms array if already build or build it
430
431 if(fContainer) return fContainer;
432 fContainer = new TObjArray(kNtasks);
433 //fContainer->SetOwner(kTRUE);
434
4226db3e 435 TH3S *h3 = NULL;
436 TObjArray *arr = NULL;
1ee39b3a 437
82b61d3c 438 // add resolution/pulls plots for dydx=ExB
439 fContainer->AddAt(arr = new TObjArray(2), kCenter);
1ee39b3a 440 arr->SetName("Center");
82b61d3c 441 if(!(h3=(TH3S*)gROOT->FindObject(Form("hCenRes%03d",fDet)))) {
442 h3 = new TH3S(
443 Form("hCenRes%03d",fDet),
444 Form(" Det[%d] Col[%d] Row[%d];t [bin];y [pw];#Delta y[cm]", fDet, fCol, fRow),
445 AliTRDseedV1::kNtb, -.5, AliTRDseedV1::kNtb-0.5, // x
446 51, -.51, .51, // y
447 60, -fDyRange, fDyRange); // dy
448 } h3->Reset();
449 arr->AddAt(h3, 0);
450 // add Pull plot for each layer
451 if(!(h3=(TH3S*)gROOT->FindObject(Form("hCenPull%03d", fDet)))){
452 h3 = new TH3S(
453 Form("hCenPull%03d", fDet),
454 Form(" Det[%d] Col[%d] Row[%d];t [bin];y [pw];#Delta y[cm]/#sigma_{y}", fDet, fCol, fRow),
455 AliTRDseedV1::kNtb, -0.5, AliTRDseedV1::kNtb-0.5, // x
456 51, -.51, .51, // y
457 60, -4., 4.); // dy/sy
458 } h3->Reset();
459 arr->AddAt(h3, 1);
460
461 if(!(h3 = (TH3S*)gROOT->FindObject(Form("Charge%03d", fDet)))){
462 h3 = new TH3S(Form("Charge%03d", fDet),
463 "dy=f(q);log(q) [a.u.];#Delta y[cm];#Delta y/#sigma_{y}",
464 50, 2.2, 7.5, 60, -fDyRange, fDyRange, 60, -4., 4.);
1ee39b3a 465 }
466 fContainer->AddAt(h3, kQRes);
467
5935a6da 468 fContainer->AddAt(arr = new TObjArray(AliTRDseedV1::kNtb), kSigm);
1ee39b3a 469 arr->SetName("Resolution");
82b61d3c 470 for(Int_t it=0; it<AliTRDseedV1::kNtb; it++){
471 if(!(h3=(TH3S*)gROOT->FindObject(Form("hr%03d_t%02d", fDet, it)))){
1ee39b3a 472 h3 = new TH3S(
7bf75b68 473 Form("hr%03d_t%02d", fDet, it),
82b61d3c 474 Form(" Det[%d] t_{drift}(%2d)[bin];z [mm];tg#phi;#Delta y[cm]", fDet, it),
1ee39b3a 475 kND, 0., 2.5, // z
476 35, -.35, .35, // tgp
563d1b38 477 60, -fDyRange, fDyRange); // dy
1ee39b3a 478 }
82b61d3c 479 arr->AddAt(h3, it);
1ee39b3a 480 }
481
5935a6da 482 fContainer->AddAt(arr = new TObjArray(AliTRDseedV1::kNtb), kMean);
1ee39b3a 483 arr->SetName("Systematics");
82b61d3c 484 for(Int_t it=0; it<AliTRDseedV1::kNtb; it++){
485 if(!(h3=(TH3S*)gROOT->FindObject(Form("hs%03d_t%02d", fDet, it)))){
1ee39b3a 486 h3 = new TH3S(
82b61d3c 487 Form("hs%03d_t%02d", fDet, it),
488 Form(" Det[%d] t_{drift}(%2d)[bin];z [mm];tg#phi - h*tg(#theta);#Delta y[cm]", fDet, it),
1ee39b3a 489 kND, 0., 2.5, // z
490 35, -.35, .35, // tgp-h tgt
563d1b38 491 60, -fDyRange, fDyRange); // dy
1ee39b3a 492 }
82b61d3c 493 arr->AddAt(h3, it);
1ee39b3a 494 }
495
496 return fContainer;
497}
498
499//_______________________________________________________
f8f46e4d 500void AliTRDclusterResolution::UserExec(Option_t *)
1ee39b3a 501{
502// Fill container histograms
503
e3147647 504
5935a6da 505 fInfo = dynamic_cast<TObjArray *>(GetInputData(1));
3c5f2bf7 506 AliDebug(2, Form("Clusters[%d]", fInfo->GetEntriesFast()));
e3147647 507 if(!IsCalibrated()){
508 LoadCalibration();
509 if(!IsCalibrated()){
510 AliWarning("Loading the calibration settings failed. Check OCDB access.");
801d4d50 511 return;
512 }
513 }
1ee39b3a 514
515 Int_t det, t;
516 Float_t x, y, z, q, dy, dydx, dzdx, cov[3], covcl[3];
4226db3e 517 TH3S *h3 = NULL;
1ee39b3a 518
519 // define limits around ExB for which x contribution is negligible
520 const Float_t kDtgPhi = 3.5e-2; //(+- 2 deg)
521
522 TObjArray *arr0 = (TObjArray*)fContainer->At(kCenter);
523 TObjArray *arr1 = (TObjArray*)fContainer->At(kSigm);
524 TObjArray *arr2 = (TObjArray*)fContainer->At(kMean);
525
4226db3e 526 const AliTRDclusterInfo *cli = NULL;
1ee39b3a 527 TIterator *iter=fInfo->MakeIterator();
528 while((cli=dynamic_cast<AliTRDclusterInfo*>((*iter)()))){
529 cli->GetCluster(det, x, y, z, q, t, covcl);
5935a6da 530
801d4d50 531 // select cluster according to detector region if specified
1ee39b3a 532 if(fDet>=0 && fDet!=det) continue;
801d4d50 533 if(fCol>=0 && fRow>=0){
534 Int_t c,r;
535 cli->GetCenterPad(c, r);
536 if(TMath::Abs(fCol-c) > 5) continue;
537 if(TMath::Abs(fRow-r) > 2) continue;
538 }
801d4d50 539 dy = cli->GetResolution();
540 AliDebug(4, Form("det[%d] tb[%2d] q[%4.0f Log[%6.4f]] dy[%7.2f][um] ypull[%5.2f]", det, t, q, TMath::Log(q), 1.e4*dy, dy/TMath::Sqrt(covcl[0])));
1ee39b3a 541
1ee39b3a 542 cli->GetGlobalPosition(y, z, dydx, dzdx, &cov[0]);
543
544 // resolution as a function of cluster charge
545 // only for phi equal exB
546 if(TMath::Abs(dydx-fExB) < kDtgPhi){
547 h3 = (TH3S*)fContainer->At(kQRes);
548 h3->Fill(TMath::Log(q), dy, dy/TMath::Sqrt(covcl[0]));
1ee39b3a 549 }
550
551 // do not use problematic clusters in resolution analysis
552 // TODO define limits as calibration aware (gain) !!
e3147647 553 if(q<20.*fGain || q>250.*fGain) continue;
1ee39b3a 554
5935a6da 555 //x = (t+.5)*fgkTimeBinLength; // conservative approach !!
1ee39b3a 556
557 // resolution as a function of y displacement from pad center
558 // only for phi equal exB
5935a6da 559 if(TMath::Abs(dydx-fExB) < kDtgPhi){
82b61d3c 560 h3 = (TH3S*)arr0->At(0);
5935a6da 561 h3->Fill(t, cli->GetYDisplacement(), dy);
82b61d3c 562 h3 = (TH3S*)arr0->At(1);
5935a6da 563 h3->Fill(t, cli->GetYDisplacement(), dy/TMath::Sqrt(covcl[0]));
1ee39b3a 564 }
565
5935a6da 566 Int_t it(((TH3S*)arr0->At(0))->GetXaxis()->FindBin(t));
1ee39b3a 567
568 // fill histo for resolution (sigma)
5935a6da 569 ((TH3S*)arr1->At(it-1))->Fill(10.*cli->GetAnisochronity(), dydx, dy);
1ee39b3a 570
571 // fill histo for systematic (mean)
5935a6da 572 ((TH3S*)arr2->At(it-1))->Fill(10.*cli->GetAnisochronity(), dydx-cli->GetTilt()*dzdx, dy);
1ee39b3a 573 }
1ee39b3a 574}
575
576
577//_______________________________________________________
578Bool_t AliTRDclusterResolution::PostProcess()
579{
64d57299 580// Steer processing of various cluster resolution dependences :
581//
582// - process resolution dependency cluster charge
583// if(HasProcess(kQRes)) ProcessCharge();
584// - process resolution dependency on y displacement
585// if(HasProcess(kCenter)) ProcessCenterPad();
586// - process resolution dependency on drift legth and drift cell width
587// if(HasProcess(kSigm)) ProcessSigma();
588// - process systematic shift on drift legth and drift cell width
589// if(HasProcess(kMean)) ProcessMean();
590
1ee39b3a 591 if(!fContainer) return kFALSE;
e3147647 592 if(!IsCalibrated()){
593 AliWarning("Not calibrated.");
594 return kFALSE;
595 }
4226db3e 596 TObjArray *arr = NULL;
597 TTree *t=NULL;
1ee39b3a 598 if(!fResults){
4226db3e 599 TGraphErrors *g = NULL;
1ee39b3a 600 fResults = new TObjArray(kNtasks);
601 fResults->SetOwner();
602 fResults->AddAt(arr = new TObjArray(3), kQRes);
603 arr->SetOwner();
604 arr->AddAt(g = new TGraphErrors(), 0);
605 g->SetLineColor(kBlue); g->SetMarkerColor(kBlue);
606 g->SetMarkerStyle(7);
607 arr->AddAt(g = new TGraphErrors(), 1);
608 g->SetLineColor(kRed); g->SetMarkerColor(kRed);
609 g->SetMarkerStyle(23);
610 arr->AddAt(g = new TGraphErrors(), 2);
611 g->SetLineColor(kGreen); g->SetMarkerColor(kGreen);
612 g->SetMarkerStyle(7);
613
614 // pad center dependence
615 fResults->AddAt(arr = new TObjArray(AliTRDgeometry::kNlayer+1), kCenter);
616 arr->SetOwner();
617 arr->AddAt(
618 t = new TTree("cent", "dy=f(y,x,ly)"), 0);
619 t->Branch("ly", &fLy, "ly/B");
5935a6da 620 t->Branch("t", &fT, "t/F");
1ee39b3a 621 t->Branch("y", &fY, "y/F");
622 t->Branch("m", &fR[0], "m[2]/F");
623 t->Branch("s", &fR[2], "s[2]/F");
624 t->Branch("pm", &fP[0], "pm[2]/F");
625 t->Branch("ps", &fP[2], "ps[2]/F");
626 for(Int_t il=1; il<=AliTRDgeometry::kNlayer; il++){
627 arr->AddAt(g = new TGraphErrors(), il);
628 g->SetLineColor(il); g->SetLineStyle(il);
629 g->SetMarkerColor(il);g->SetMarkerStyle(4);
630 }
631
632
633 fResults->AddAt(t = new TTree("sigm", "dy=f(dw,x,dydx)"), kSigm);
5935a6da 634 t->Branch("t", &fT, "t/F");
635 t->Branch("x", &fX, "x/F");
1ee39b3a 636 t->Branch("z", &fZ, "z/F");
637 t->Branch("sx", &fR[0], "sx[2]/F");
638 t->Branch("sy", &fR[2], "sy[2]/F");
639
640
641 fResults->AddAt(t = new TTree("mean", "dy=f(dw,x,dydx - h dzdx)"), kMean);
5935a6da 642 t->Branch("t", &fT, "t/F");
643 t->Branch("x", &fX, "x/F");
1ee39b3a 644 t->Branch("z", &fZ, "z/F");
645 t->Branch("dx", &fR[0], "dx[2]/F");
646 t->Branch("dy", &fR[2], "dy[2]/F");
647 } else {
4226db3e 648 TObject *o = NULL;
1ee39b3a 649 TIterator *iter=fResults->MakeIterator();
650 while((o=((*iter)()))) o->Clear(); // maybe it is wrong but we should never reach this point
651 }
652
1ee39b3a 653 // process resolution dependency on charge
654 if(HasProcess(kQRes)) ProcessCharge();
655
656 // process resolution dependency on y displacement
657 if(HasProcess(kCenter)) ProcessCenterPad();
658
659 // process resolution dependency on drift legth and drift cell width
660 if(HasProcess(kSigm)) ProcessSigma();
661
662 // process systematic shift on drift legth and drift cell width
663 if(HasProcess(kMean)) ProcessMean();
664
665 return kTRUE;
666}
667
668//_______________________________________________________
e3147647 669Bool_t AliTRDclusterResolution::LoadCalibration()
1ee39b3a 670{
801d4d50 671// Retrieve calibration parameters from OCDB, drift velocity and t0 for the detector region specified by
672// a previous call to AliTRDclusterResolution::SetCalibrationRegion().
673
56f313bd 674 AliCDBManager *cdb = AliCDBManager::Instance(); // check access OCDB
1ee39b3a 675 if(cdb->GetRun() < 0){
676 AliError("OCDB manager not properly initialized");
677 return kFALSE;
678 }
1ee39b3a 679 // check magnetic field
56f313bd 680 if(!TGeoGlobalMagField::Instance() || !TGeoGlobalMagField::Instance()->IsLocked()){
681 AliError("Magnetic field not available.");
801d4d50 682 return kFALSE;
1ee39b3a 683 }
684
801d4d50 685 // check pad for detector
686 if(fCol>=0 && fRow>=0){
687 AliTRDgeometry geo;
688 AliTRDpadPlane *pp(geo.GetPadPlane(fDet));
689 if(fCol>=pp->GetNcols() ||
690 fRow>=pp->GetNrows()){
691 AliWarning(Form("Pad coordinates col[%d] or row[%d] incorrect for det[%d].\nLimits are max col[%d] max row[%d]. Reset to default", fCol, fRow, fDet, pp->GetNcols(), pp->GetNrows()));
692 fCol = -1; fRow=-1;
693 }
694 }
1ee39b3a 695
696 AliTRDcalibDB *fCalibration = AliTRDcalibDB::Instance();
801d4d50 697 AliTRDCalROC *fCalVdriftROC(fCalibration->GetVdriftROC(fDet>=0?fDet:0))
698 ,*fCalT0ROC(fCalibration->GetT0ROC(fDet>=0?fDet:0));
1ee39b3a 699 const AliTRDCalDet *fCalVdriftDet = fCalibration->GetVdriftDet();
5935a6da 700 const AliTRDCalDet *fCalT0Det = fCalibration->GetT0Det();
1ee39b3a 701
56f313bd 702 if(IsUsingCalibParam(kVdrift)){
703 fVdrift = fCalVdriftDet->GetValue(fDet>=0?fDet:0);
704 if(fCol>=0 && fRow>=0) fVdrift*= fCalVdriftROC->GetValue(fCol, fRow);
705 }
5935a6da 706 fExB = AliTRDCommonParam::Instance()->GetOmegaTau(fVdrift);
56f313bd 707 if(IsUsingCalibParam(kT0)){
708 fT0 = fCalT0Det->GetValue(fDet>=0?fDet:0);
709 if(fCol>=0 && fRow>=0) fT0 *= fCalT0ROC->GetValue(fCol, fRow);
710 }
711 if(IsUsingCalibParam(kGain)) fGain = (fCol>=0 && fRow>=0)?fCalibration-> GetGainFactor(fDet, fCol, fRow):fCalibration-> GetGainFactorAverage(fDet);
712
e3147647 713 SetBit(kCalibrated);
5935a6da 714
56f313bd 715 AliInfo(Form("Calibrate for Det[%3d] Col[%3d] Row[%2d] : \n t0[%5.3f] vd[%5.3f] gain[%5.3f] ExB[%f]", fDet, fCol, fRow, fT0, fVdrift, fGain, fExB));
5935a6da 716
1ee39b3a 717 return kTRUE;
718}
719
801d4d50 720//_______________________________________________________
721void AliTRDclusterResolution::SetCalibrationRegion(Int_t det, Int_t col, Int_t row)
722{
723// Set calibration region in terms of detector and pad.
724// By default detector 0 mean values are considered.
725
726 if(det>=0 && det<AliTRDgeometry::kNdet){
727 fDet = det;
728 if(col>=0 && row>=0){
729 fCol = col;
730 fRow = row;
731 }
732 return;
733 }
734 AliError(Form("Detector index outside range [0 %d].", AliTRDgeometry::kNdet));
735}
736
1ee39b3a 737//_______________________________________________________
738void AliTRDclusterResolution::SetVisual()
739{
740 if(fCanvas) return;
741 fCanvas = new TCanvas("clResCanvas", "Cluster Resolution Visualization", 10, 10, 600, 600);
742}
743
744//_______________________________________________________
745void AliTRDclusterResolution::ProcessCharge()
746{
747// Resolution as a function of cluster charge.
748//
749// As described in the function ProcessCenterPad() the error parameterization for clusters for phi = a_L can be
750// written as:
751// BEGIN_LATEX
752// #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
753// END_LATEX
754// with the contribution in case of B=0 given by:
755// BEGIN_LATEX
756// #sigma_{y}|_{B=0} = #sigma_{diff}*Gauss(0, s_{ly}) + #delta_{#sigma}(q)
757// END_LATEX
758// which further can be simplified to:
759// BEGIN_LATEX
760// <#sigma_{y}|_{B=0}>(q) = <#sigma_{y}> + #delta_{#sigma}(q)
761// <#sigma_{y}> = #int{f(q)#sigma_{y}dq}
762// END_LATEX
763// The results for s_y and f(q) are displayed below:
764//Begin_Html
765//<img src="TRD/clusterQerror.gif">
766//End_Html
767// The function has to extended to accomodate gain calibration scalling and errors.
768//
769// Author
770// Alexandru Bercuci <A.Bercuci@gsi.de>
771
5935a6da 772 TH3S *h3(NULL);
773 if(!(h3 = (TH3S*)fContainer->At(kQRes))) {
1ee39b3a 774 AliWarning("Missing dy=f(Q) histo");
775 return;
776 }
777 TF1 f("f", "gaus", -.5, .5);
5935a6da 778 TAxis *ax(NULL);
779 TH1 *h1(NULL);
1ee39b3a 780
781 // compute mean error on x
782 Double_t s2x = 0.;
5935a6da 783 for(Int_t ix=5; ix<AliTRDseedV1::kNtb; ix++){
1ee39b3a 784 // retrieve error on the drift length
785 s2x += AliTRDcluster::GetSX(ix);
786 }
5935a6da 787 s2x /= (AliTRDseedV1::kNtb-5); s2x *= s2x;
76d976d2 788 //Double_t exb2 = fExB*fExB;
1ee39b3a 789
790 TObjArray *arr = (TObjArray*)fResults->At(kQRes);
791 TGraphErrors *gqm = (TGraphErrors*)arr->At(0);
792 TGraphErrors *gqs = (TGraphErrors*)arr->At(1);
793 TGraphErrors *gqp = (TGraphErrors*)arr->At(2);
794 Double_t q, n = 0., entries;
5935a6da 795 ax = h3->GetXaxis();
1ee39b3a 796 for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
797 q = TMath::Exp(ax->GetBinCenter(ix));
5935a6da 798 ax->SetRange(ix, ix);
799 h1 = h3->Project3D("y");
1ee39b3a 800 entries = h1->GetEntries();
5935a6da 801 if(entries < 150) continue;
1ee39b3a 802 h1->Fit(&f, "Q");
803
804 // Fill sy^2 = f(q)
805 Int_t ip = gqm->GetN();
806 gqm->SetPoint(ip, q, 1.e4*f.GetParameter(1));
807 gqm->SetPointError(ip, 0., 1.e4*f.GetParError(1));
808
809 // correct sigma for ExB effect
5935a6da 810 gqs->SetPoint(ip, q, 1.e4*f.GetParameter(2)/**f.GetParameter(2)-exb2*s2x)*/);
811 gqs->SetPointError(ip, 0., 1.e4*f.GetParError(2)/**f.GetParameter(2)*/);
1ee39b3a 812
813 // save probability
814 n += entries;
815 gqp->SetPoint(ip, q, entries);
816 gqp->SetPointError(ip, 0., 0./*TMath::Sqrt(entries)*/);
817 }
818
819 // normalize probability and get mean sy
820 Double_t sm = 0., sy;
821 for(Int_t ip=gqp->GetN(); ip--;){
822 gqp->GetPoint(ip, q, entries);
823 entries/=n;
5935a6da 824 gqp->SetPoint(ip, q, 1.e4*entries);
1ee39b3a 825 gqs->GetPoint(ip, q, sy);
826 sm += entries*sy;
827 }
828
829 // error parametrization s(q) = <sy> + b(1/q-1/q0)
830 TF1 fq("fq", "[0] + [1]/x", 20., 250.);
831 gqs->Fit(&fq/*, "W"*/);
832 printf("sm=%f [0]=%f [1]=%f\n", 1.e-4*sm, fq.GetParameter(0), fq.GetParameter(1));
833 printf(" const Float_t sq0inv = %f; // [1/q0]\n", (sm-fq.GetParameter(0))/fq.GetParameter(1));
834 printf(" const Float_t sqb = %f; // [cm]\n", 1.e-4*fq.GetParameter(1));
835}
836
837//_______________________________________________________
838void AliTRDclusterResolution::ProcessCenterPad()
839{
840// Resolution as a function of y displacement from pad center and drift length.
841//
842// Since the error parameterization of cluster r-phi position can be written as (see AliTRDcluster::SetSigmaY2()):
843// BEGIN_LATEX
844// #sigma_{y}^{2} = (#sigma_{diff}*Gauss(0, s_{ly}) + #delta_{#sigma}(q))^{2} + tg^{2}(#alpha_{L})*#sigma_{x}^{2} + tg^{2}(#phi-#alpha_{L})*#sigma_{x}^{2}+[tg(#phi-#alpha_{L})*tg(#alpha_{L})*x]^{2}/12
845// END_LATEX
846// one can see that for phi = a_L one gets the following expression:
847// BEGIN_LATEX
848// #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
849// END_LATEX
850// where we have explicitely marked the remaining term in case of absence of magnetic field. Thus one can use the
851// previous equation to estimate s_y for B=0 and than by comparing in magnetic field conditions one can get the s_x.
852// This is a simplified method to determine the error parameterization for s_x and s_y as compared to the one
853// implemented in ProcessSigma(). For more details on cluster error parameterization please see also
854// AliTRDcluster::SetSigmaY2()
855//
856// The representation of dy=f(y_cen, x_drift| layer) can be also used to estimate the systematic shift in the r-phi
857// coordinate resulting from imperfection in the cluster shape parameterization. From the expresion of the shift derived
858// in ProcessMean() with phi=exb one gets:
859// BEGIN_LATEX
860// <#Delta y>= <#delta x> * (tg(#alpha_{L})-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
861// <#Delta y>(y_{cen})= -h*<#delta x>(x_{drift}, q_{cl}) * dz/dx + #delta y(y_{cen}, ...)
862// END_LATEX
863// where all dependences are made explicit. This last expression can be used in two ways:
864// - by average on the dz/dx we can determine directly dy (the method implemented here)
865// - by plotting as a function of dzdx one can determine both dx and dy components in an independent method.
866//Begin_Html
867//<img src="TRD/clusterYcorr.gif">
868//End_Html
869// Author
870// Alexandru Bercuci <A.Bercuci@gsi.de>
871
872 TObjArray *arr = (TObjArray*)fContainer->At(kCenter);
873 if(!arr) {
874 AliWarning("Missing dy=f(y | x, ly) container");
875 return;
876 }
877 Double_t exb2 = fExB*fExB;
878 Float_t s[AliTRDgeometry::kNlayer];
879 TF1 f("f", "gaus", -.5, .5);
880 TF1 fp("fp", "gaus", -3.5, 3.5);
881
4226db3e 882 TH1D *h1 = NULL; TH2F *h2 = NULL; TH3S *h3r=NULL, *h3p=NULL;
1ee39b3a 883 TObjArray *arrRes = (TObjArray*)fResults->At(kCenter);
884 TTree *t = (TTree*)arrRes->At(0);
4226db3e 885 TGraphErrors *gs = NULL;
886 TAxis *ax = NULL;
1ee39b3a 887
5935a6da 888 AliDebug(1, Form("Calibrate for Det[%3d] t0[%5.3f] vd[%5.3f]", fDet, fT0, fVdrift));
889
1ee39b3a 890 const Int_t nl = AliTRDgeometry::kNlayer;
5935a6da 891 printf(" const Float_t lSy[%d][%d] = {\n {", nl, AliTRDseedV1::kNtb);
1ee39b3a 892 for(Int_t il=0; il<nl; il++){
893 if(!(h3r = (TH3S*)arr->At(il))) continue;
894 if(!(h3p = (TH3S*)arr->At(nl+il))) continue;
895 gs = (TGraphErrors*)arrRes->At(il+1);
896 fLy = il;
1ee39b3a 897 for(Int_t ix=1; ix<=h3r->GetXaxis()->GetNbins(); ix++){
898 ax = h3r->GetXaxis(); ax->SetRange(ix, ix);
899 ax = h3p->GetXaxis(); ax->SetRange(ix, ix);
5935a6da 900 fT = ax->GetBinCenter(ix);
1ee39b3a 901 for(Int_t iy=1; iy<=h3r->GetYaxis()->GetNbins(); iy++){
902 ax = h3r->GetYaxis(); ax->SetRange(iy, iy);
903 ax = h3p->GetYaxis(); ax->SetRange(iy, iy);
904 fY = ax->GetBinCenter(iy);
1ee39b3a 905 // finish navigation in the HnSparse
906
907 h1 = (TH1D*)h3r->Project3D("z");
908 Int_t entries = (Int_t)h1->Integral();
909 if(entries < 50) continue;
910 //Adjust(&f, h1);
911 h1->Fit(&f, "QN");
912
913 // Fill sy,my=f(y_w,x,ly)
914 fR[0] = f.GetParameter(1); fR[1] = f.GetParError(1);
915 fR[2] = f.GetParameter(2); fR[3] = f.GetParError(2);
916
917 h1 = (TH1D*)h3p->Project3D("z");
918 h1->Fit(&fp, "QN");
919 fP[0] = fp.GetParameter(1); fP[1] = fp.GetParError(1);
920 fP[2] = fp.GetParameter(2); fP[3] = fp.GetParError(2);
921
ca50a37e 922 AliDebug(4, Form("ly[%d] tb[%2d] y[%+5.2f] m[%5.3f] s[%5.3f] pm[%5.3f] ps[%5.3f]", fLy, (Int_t)fT, fY, fR[0], fR[2], fP[0], fP[2]));
1ee39b3a 923 t->Fill();
1ee39b3a 924 }
925 }
5935a6da 926 t->Draw(Form("y:t>>h(%d, -0.5, %f, 51, -.51, .51)", AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5),
1ee39b3a 927 Form("s[0]*(ly==%d&&abs(m[0])<1.e-1)", fLy),
928 "goff");
929 h2=(TH2F*)gROOT->FindObject("h");
930 f.FixParameter(1, 0.);
931 Int_t n = h2->GetXaxis()->GetNbins(), nn(0); s[il]=0.;
932 printf(" {");
933 for(Int_t ix=1; ix<=n; ix++){
934 ax = h2->GetXaxis();
5935a6da 935 fT = ax->GetBinCenter(ix);
1ee39b3a 936 h1 = h2->ProjectionY("hCenPy", ix, ix);
937 //if((Int_t)h1->Integral() < 1.e-10) continue;
938
939 // Apply lorentz angle correction
940 // retrieve error on the drift length
941 Double_t s2x = AliTRDcluster::GetSX(ix-1); s2x *= s2x;
942 Int_t nnn = 0;
943 for(Int_t iy=1; iy<=h1->GetNbinsX(); iy++){
944 Double_t s2 = h1->GetBinContent(iy); s2*= s2;
945 // sigma square corrected for Lorentz angle
946 // s2 = s2_y(y_w,x)+exb2*s2_x
947 Double_t sy = TMath::Sqrt(TMath::Max(s2 - exb2*s2x, Double_t(0.)));
948 if(sy<1.e-20) continue;
949 h1->SetBinContent(iy, sy); nnn++;
5935a6da 950 AliDebug(4, Form("s[%6.2f] sx[%6.2f] sy[%6.2f]\n",
1ee39b3a 951 1.e4*TMath::Sqrt(s2), 1.e4*TMath::Abs(fExB*AliTRDcluster::GetSX(ix-1)),
5935a6da 952 1.e4*h1->GetBinContent(iy)));
1ee39b3a 953 }
954 // do fit only if enough data
955 Double_t sPRF = 0.;
956 if(nnn>5){
957 h1->Fit(&f, "QN");
5935a6da 958 sPRF = f.GetParameter(2); nn++;
1ee39b3a 959 }
960 s[il]+=sPRF;
961 printf("%6.4f,%s", sPRF, ix%6?" ":"\n ");
962 Int_t jx = gs->GetN();
5935a6da 963 gs->SetPoint(jx, fT, 1.e4*sPRF);
1ee39b3a 964 gs->SetPointError(jx, 0., 0./*f.GetParError(0)*/);
965 }
966 printf("\b},\n");
967 s[il]/=nn;
968
5935a6da 969 f.ReleaseParameter(1);
1ee39b3a 970
971
972 if(!fCanvas) continue;
973 h2->Draw("lego2fb");
974 fCanvas->Modified(); fCanvas->Update();
975 if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessCenter_ly[%d].gif", fLy));
976 else gSystem->Sleep(100);
977 }
978 printf(" };\n");
979 printf(" const Float_t lPRF[] = {"
980 "%5.3f, %5.3f, %5.3f, %5.3f, %5.3f, %5.3f};\n",
981 s[0], s[1], s[2], s[3], s[4], s[5]);
982}
983
984//_______________________________________________________
985void AliTRDclusterResolution::ProcessSigma()
986{
987// As the r-phi coordinate is the only one which is measured by the TRD detector we have to rely on it to
988// estimate both the radial (x) and r-phi (y) errors. This method is based on the following assumptions.
989// The measured error in the y direction is the sum of the intrinsic contribution of the r-phi measurement
990// with the contribution of the radial measurement - because x is not a parameter of Alice track model (Kalman).
991// BEGIN_LATEX
992// #sigma^{2}|_{y} = #sigma^{2}_{y*} + #sigma^{2}_{x*}
993// END_LATEX
994// In the general case
995// BEGIN_LATEX
996// #sigma^{2}_{y*} = #sigma^{2}_{y} + tg^{2}(#alpha_{L})#sigma^{2}_{x_{drift}}
997// #sigma^{2}_{x*} = tg^{2}(#phi - #alpha_{L})*(#sigma^{2}_{x_{drift}} + #sigma^{2}_{x_{0}} + tg^{2}(#alpha_{L})*x^{2}/12)
998// END_LATEX
999// where we have explicitely show the lorentz angle correction on y and the projection of radial component on the y
1000// direction through the track angle in the bending plane (phi). Also we have shown that the radial component in the
1001// last equation has twp terms, the drift and the misalignment (x_0). For ideal geometry or known misalignment one
1002// can solve the equation
1003// BEGIN_LATEX
1004// #sigma^{2}|_{y} = tg^{2}(#phi - #alpha_{L})*(#sigma^{2}_{x} + tg^{2}(#alpha_{L})*x^{2}/12)+ [#sigma^{2}_{y} + tg^{2}(#alpha_{L})#sigma^{2}_{x}]
1005// END_LATEX
1006// by fitting a straight line:
1007// BEGIN_LATEX
1008// #sigma^{2}|_{y} = a(x_{cl}, z_{cl}) * tg^{2}(#phi - #alpha_{L}) + b(x_{cl}, z_{cl})
1009// END_LATEX
1010// the error parameterization will be given by:
1011// BEGIN_LATEX
1012// #sigma_{x} (x_{cl}, z_{cl}) = #sqrt{a(x_{cl}, z_{cl}) - tg^{2}(#alpha_{L})*x^{2}/12}
1013// #sigma_{y} (x_{cl}, z_{cl}) = #sqrt{b(x_{cl}, z_{cl}) - #sigma^{2}_{x} (x_{cl}, z_{cl}) * tg^{2}(#alpha_{L})}
1014// END_LATEX
1015// Below there is an example of such dependency.
1016//Begin_Html
1017//<img src="TRD/clusterSigmaMethod.gif">
1018//End_Html
1019//
1020// The error parameterization obtained by this method are implemented in the functions AliTRDcluster::GetSX() and
1021// AliTRDcluster::GetSYdrift(). For an independent method to determine s_y as a function of drift length check the
1022// function ProcessCenterPad(). One has to keep in mind that while this method return the mean s_y over the distance
1023// to pad center distribution the other method returns the *STANDARD* value at center=0 (maximum). To recover the
1024// standard value one has to solve the obvious equation:
1025// BEGIN_LATEX
1026// #sigma_{y}^{STANDARD} = #frac{<#sigma_{y}>}{#int{s exp(s^{2}/#sigma) ds}}
1027// END_LATEX
1028// with "<s_y>" being the value calculated here and "sigma" the width of the s_y distribution calculated in
1029// ProcessCenterPad().
1030//
1031// Author
1032// Alexandru Bercuci <A.Bercuci@gsi.de>
1033
1034 TObjArray *arr = (TObjArray*)fContainer->At(kSigm);
1035 if(!arr){
1036 AliWarning("Missing dy=f(x_d, d_w) container");
1037 return;
1038 }
1039
1040 // init visualization
4226db3e 1041 TGraphErrors *ggs = NULL;
1042 TGraph *line = NULL;
1ee39b3a 1043 if(fCanvas){
1044 ggs = new TGraphErrors();
1045 line = new TGraph();
1046 line->SetLineColor(kRed);line->SetLineWidth(2);
1047 }
1048
1049 // init logistic support
1050 TF1 f("f", "gaus", -.5, .5);
1051 TLinearFitter gs(1,"pol1");
4226db3e 1052 TH1 *hFrame=NULL;
1053 TH1D *h1 = NULL; TH3S *h3=NULL;
1054 TAxis *ax = NULL;
5935a6da 1055 Double_t exb2 = fExB*fExB;
1ee39b3a 1056 AliTRDcluster c;
1057 TTree *t = (TTree*)fResults->At(kSigm);
5935a6da 1058 for(Int_t ix=0; ix<AliTRDseedV1::kNtb; ix++){
1ee39b3a 1059 if(!(h3=(TH3S*)arr->At(ix))) continue;
1060 c.SetPadTime(ix);
5935a6da 1061 fX = c.GetXloc(fT0, fVdrift);
1062 fT = c.GetLocalTimeBin(); // ideal
1063 printf(" pad time[%d] local[%f]\n", ix, fT);
1ee39b3a 1064 for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
1065 ax = h3->GetXaxis();
1066 ax->SetRange(iz, iz);
1067 fZ = ax->GetBinCenter(iz);
1068
1069 // reset visualization
1070 if(fCanvas){
1071 new(ggs) TGraphErrors();
1072 ggs->SetMarkerStyle(7);
1073 }
1074 gs.ClearPoints();
1075
1076 for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
1077 ax = h3->GetYaxis();
1078 ax->SetRange(ip, ip);
1079 Double_t tgl = ax->GetBinCenter(ip);
1080 // finish navigation in the HnSparse
1081
1082 //if(TMath::Abs(dydx)>0.18) continue;
1083 Double_t tgg = (tgl-fExB)/(1.+tgl*fExB);
1084 Double_t tgg2 = tgg*tgg;
1085
1086 h1 = (TH1D*)h3->Project3D("z");
1087 Int_t entries = (Int_t)h1->Integral();
1088 if(entries < 50) continue;
1089 //Adjust(&f, h1);
1090 h1->Fit(&f, "QN");
1091
1092 Double_t s2 = f.GetParameter(2)*f.GetParameter(2);
1093 Double_t s2e = 2.*f.GetParameter(2)*f.GetParError(2);
1094 // Fill sy^2 = f(tg^2(phi-a_L))
1095 gs.AddPoint(&tgg2, s2, s2e);
1096
1097 if(!ggs) continue;
1098 Int_t jp = ggs->GetN();
1099 ggs->SetPoint(jp, tgg2, s2);
1100 ggs->SetPointError(jp, 0., s2e);
1101 }
1102 // TODO here a more robust fit method has to be provided
1103 // for which lower boundaries on the parameters have to
1104 // be imposed. Unfortunately the Minuit fit does not work
1105 // for the TGraph in the case of B not 0.
1106 if(gs.Eval()) continue;
1107
5935a6da 1108 fR[0] = gs.GetParameter(1) - fX*fX*exb2/12.;
1109 AliDebug(3, Form(" s2x+x2=%f ang=%f s2x=%f", gs.GetParameter(1), fX*fX*exb2/12., fR[0]));
1ee39b3a 1110 fR[0] = TMath::Max(fR[0], Float_t(4.e-4));
1111
1112 // s^2_y = s0^2_y + tg^2(a_L) * s^2_x
1113 // s0^2_y = f(D_L)*x + s_PRF^2
1114 fR[2]= gs.GetParameter(0)-exb2*fR[0];
5935a6da 1115 AliDebug(3, Form(" s2y+s2x=%f s2y=%f", fR[0], fR[2]));
1ee39b3a 1116 fR[2] = TMath::Max(fR[2], Float_t(2.5e-5));
1117 fR[0] = TMath::Sqrt(fR[0]);
1118 fR[1] = .5*gs.GetParError(1)/fR[0];
1119 fR[2] = TMath::Sqrt(fR[2]);
1120 fR[3] = gs.GetParError(0)+exb2*exb2*gs.GetParError(1);
1121 t->Fill();
5935a6da 1122 AliDebug(2, Form("xd=%4.2f[cm] sx=%6.1f[um] sy=%5.1f[um]", fX, 1.e4*fR[0], 1.e4*fR[2]));
1ee39b3a 1123
1124 if(!fCanvas) continue;
1125 fCanvas->cd(); fCanvas->SetLogx(); //fCanvas->SetLogy();
1126 if(!hFrame){
1127 fCanvas->SetMargin(0.15, 0.01, 0.1, 0.01);
1128 hFrame=new TH1I("hFrame", "", 100, 0., .3);
1129 hFrame->SetMinimum(0.);hFrame->SetMaximum(.005);
1130 hFrame->SetXTitle("tg^{2}(#phi-#alpha_{L})");
1131 hFrame->SetYTitle("#sigma^{2}y[cm^{2}]");
1132 hFrame->GetYaxis()->SetTitleOffset(2.);
1133 hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
1134 hFrame->Draw();
1135 } else hFrame->Reset();
1136 Double_t xx = 0., dxx=.2/50;
1137 for(Int_t ip=0;ip<50;ip++){
1138 line->SetPoint(ip, xx, gs.GetParameter(0)+xx*gs.GetParameter(1));
1139 xx+=dxx;
1140 }
1141 ggs->Draw("pl"); line->Draw("l");
1142 fCanvas->Modified(); fCanvas->Update();
1143 if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessSigma_z[%5.3f]_x[%5.3f].gif", fZ, fX));
1144 else gSystem->Sleep(100);
1145 }
1146 }
1147 return;
1148}
1149
1150//_______________________________________________________
1151void AliTRDclusterResolution::ProcessMean()
1152{
1153// By this method the cluster shift in r-phi and radial directions can be estimated by comparing with the MC.
1154// The resolution of the cluster corrected for pad tilt with respect to MC in the r-phi (measuring) plane can be
1155// expressed by:
1156// BEGIN_LATEX
1157// #Delta y=w - y_{MC}(x_{cl})
1158// w = y_{cl}^{'} + h*(z_{MC}(x_{cl})-z_{cl})
1159// y_{MC}(x_{cl}) = y_{0} - dy/dx*x_{cl}
1160// z_{MC}(x_{cl}) = z_{0} - dz/dx*x_{cl}
1161// y_{cl}^{'} = y_{cl}-x_{cl}*tg(#alpha_{L})
1162// END_LATEX
1163// where x_cl is the drift length attached to a cluster, y_cl is the r-phi coordinate of the cluster measured by
1164// charge sharing on adjacent pads and y_0 and z_0 are MC reference points (as example the track references at
1165// entrance/exit of a chamber). If we suppose that both r-phi (y) and radial (x) coordinate of the clusters are
1166// affected by errors we can write
1167// BEGIN_LATEX
1168// x_{cl} = x_{cl}^{*} + #delta x
1169// y_{cl} = y_{cl}^{*} + #delta y
1170// END_LATEX
1171// where the starred components are the corrected values. Thus by definition the following quantity
1172// BEGIN_LATEX
1173// #Delta y^{*}= w^{*} - y_{MC}(x_{cl}^{*})
1174// END_LATEX
1175// has 0 average over all dependency. Using this decomposition we can write:
1176// BEGIN_LATEX
1177// <#Delta y>=<#Delta y^{*}> + <#delta x * (dy/dx-h*dz/dx) + #delta y - #delta x * tg(#alpha_{L})>
1178// END_LATEX
1179// which can be transformed to the following linear dependence:
1180// BEGIN_LATEX
1181// <#Delta y>= <#delta x> * (dy/dx-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
1182// END_LATEX
1183// if expressed as function of dy/dx-h*dz/dx. Furtheremore this expression can be plotted for various clusters
1184// i.e. we can explicitely introduce the diffusion (x_cl) and drift cell - anisochronity (z_cl) dependences. From
1185// plotting this dependence and linear fitting it with:
1186// BEGIN_LATEX
1187// <#Delta y>= a(x_{cl}, z_{cl}) * (dy/dx-h*dz/dx) + b(x_{cl}, z_{cl})
1188// END_LATEX
1189// the systematic shifts will be given by:
1190// BEGIN_LATEX
1191// #delta x (x_{cl}, z_{cl}) = a(x_{cl}, z_{cl})
1192// #delta y (x_{cl}, z_{cl}) = b(x_{cl}, z_{cl}) + a(x_{cl}, z_{cl}) * tg(#alpha_{L})
1193// END_LATEX
1194// Below there is an example of such dependency.
1195//Begin_Html
1196//<img src="TRD/clusterShiftMethod.gif">
1197//End_Html
1198//
1199// The occurance of the radial shift is due to the following conditions
1200// - the approximation of a constant drift velocity over the drift length (larger drift velocities close to
1201// cathode wire plane)
1202// - the superposition of charge tails in the amplification region (first clusters appear to be located at the
1203// anode wire)
1204// - the superposition of charge tails in the drift region (shift towards anode wire)
1205// - diffusion effects which convolute with the TRF thus enlarging it
1206// - approximate knowledge of the TRF (approximate measuring in test beam conditions)
1207//
1208// The occurance of the r-phi shift is due to the following conditions
1209// - approximate model for cluster shape (LUT)
1210// - rounding-up problems
1211//
1212// The numerical results for ideal simulations for the radial and r-phi shifts are displayed below and used
1213// for the cluster reconstruction (see the functions AliTRDcluster::GetXcorr() and AliTRDcluster::GetYcorr()).
1214//Begin_Html
1215//<img src="TRD/clusterShiftX.gif">
1216//<img src="TRD/clusterShiftY.gif">
1217//End_Html
1218// More details can be found in the presentation given during the TRD
1219// software meeting at the end of 2008 and beginning of year 2009, published on indico.cern.ch.
1220//
1221// Author
1222// Alexandru Bercuci <A.Bercuci@gsi.de>
1223
1224
1225
1226 TObjArray *arr = (TObjArray*)fContainer->At(kMean);
1227 if(!arr){
1228 AliWarning("Missing dy=f(x_d, d_w) container");
1229 return;
1230 }
1231
1232 // init logistic support
1233 TF1 f("f", "gaus", -.5, .5);
1234 TF1 line("l", "[0]+[1]*x", -.15, .15);
1235 TGraphErrors *gm = new TGraphErrors();
4226db3e 1236 TH1 *hFrame=NULL;
1237 TH1D *h1 = NULL; TH3S *h3 =NULL;
1238 TAxis *ax = NULL;
5935a6da 1239
1240 AliDebug(1, Form("Calibrate for Det[%3d] t0[%5.3f] vd[%5.3f]", fDet, fT0, fVdrift));
1ee39b3a 1241
1242 AliTRDcluster c;
1243 TTree *t = (TTree*)fResults->At(kMean);
5935a6da 1244 for(Int_t ix=0; ix<AliTRDseedV1::kNtb; ix++){
1ee39b3a 1245 if(!(h3=(TH3S*)arr->At(ix))) continue;
1246 c.SetPadTime(ix);
5935a6da 1247 fX = c.GetXloc(fT0, fVdrift);
1248 fT = c.GetLocalTimeBin();
1ee39b3a 1249 for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
1250 ax = h3->GetXaxis();
1251 ax->SetRange(iz, iz);
1252 fZ = ax->GetBinCenter(iz);
1253
1254 // reset fitter
1255 new(gm) TGraphErrors();
1256 gm->SetMarkerStyle(7);
1257
1258 for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
1259 ax = h3->GetYaxis();
1260 ax->SetRange(ip, ip);
1261 Double_t tgl = ax->GetBinCenter(ip);
1262 // finish navigation in the HnSparse
1263
1264 h1 = (TH1D*)h3->Project3D("z");
1265 Int_t entries = (Int_t)h1->Integral();
b9ddd472 1266 if(entries < 50) continue;
1ee39b3a 1267 //Adjust(&f, h1);
1268 h1->Fit(&f, "QN");
1269
1270 // Fill <Dy> = f(dydx - h*dzdx)
1271 Int_t jp = gm->GetN();
1272 gm->SetPoint(jp, tgl, f.GetParameter(1));
1273 gm->SetPointError(jp, 0., f.GetParError(1));
1274 }
5935a6da 1275 if(gm->GetN()<10) continue;
1ee39b3a 1276
1277 gm->Fit(&line, "QN");
1278 fR[0] = line.GetParameter(1); // dx
1279 fR[1] = line.GetParError(1);
1280 fR[2] = line.GetParameter(0) + fExB*fR[0]; // xs = dy - tg(a_L)*dx
1281 t->Fill();
5935a6da 1282 AliDebug(2, Form("tb[%02d] xd=%4.2f[cm] dx=%6.2f[um] dy=%6.2f[um]", ix, fX, 1.e4*fR[0], 1.e4*fR[2]));
1ee39b3a 1283 if(!fCanvas) continue;
5935a6da 1284
1ee39b3a 1285 fCanvas->cd();
1286 if(!hFrame){
1287 fCanvas->SetMargin(0.1, 0.02, 0.1, 0.01);
1288 hFrame=new TH1I("hFrame", "", 100, -.3, .3);
1289 hFrame->SetMinimum(-.1);hFrame->SetMaximum(.1);
1290 hFrame->SetXTitle("tg#phi-htg#theta");
1291 hFrame->SetYTitle("#Delta y[cm]");
1292 hFrame->GetYaxis()->SetTitleOffset(1.5);
1293 hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
1294 hFrame->Draw();
1295 } else hFrame->Reset();
1296 gm->Draw("pl"); line.Draw("same");
1297 fCanvas->Modified(); fCanvas->Update();
5935a6da 1298 if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessMean_Z[%5.3f]_TB[%02d].gif", fZ, ix));
1ee39b3a 1299 else gSystem->Sleep(100);
1300 }
1301 }
1302}