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: AliTRDclusterResolution.cxx */
18 ///////////////////////////////////////////////////////////////////////////////
20 // TRD cluster error parameterization //
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//
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.
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
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.
44 // The functionalities implemented in this class are based on the storage
45 // class AliTRDclusterInfo.
50 // The method to disentangle s_y and s_x is based on the relation (see also fig.)
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}})
56 // #sigma^{2}_{x_{c}} #approx 0
58 // we suppose the chamber is well calibrated for t_{0} and aligned in
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.
67 // We estimate this effects by the relations:
69 // #mu_{y} = tg(#alpha_{L})*#Delta x_{d}(...) + tg(#phi-#alpha_{L})*(#Delta x_{c}(...) + #Delta x_{d}(...))
73 // #Delta x_{d}(...) = (<v_{d}> + #delta v_{d}(x_{d}, d)) * (t + t^{*}(Q))
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).
78 // For estimating the contribution from asymmetry of TRF the following
79 // parameterization is being used
81 // t^{*}(Q) = #delta_{0} * #frac{Q_{t+1} - Q_{t-1}}{Q_{t-1} + Q_{t} + Q_{t+1}}
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 :
89 // <#Delta y> = a + b * sin(c*y_{pw})
94 // Parameterization against total charge
96 // Obtained for B=0T at phi=0. All other effects integrated out.
98 // #sigma^{2}_{y}(Q) = #sigma^{2}_{y}(...) + b(#frac{1}{Q} - #frac{1}{Q_{0}})
100 // For B diff 0T the error of the average ExB correction error has to be subtracted !!
102 // Parameterization Sx
104 // The parameterization of the error in the x direction can be written as
106 // #sigma_{x} = #sigma_{x}^{||} + #sigma_{x}^{#perp}
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
118 // and perpendicular to the track is pol2 with the parameters
120 // Par_0 = 0.190676 +/- 0.41785
121 // Par_1 = -3.9269 +/- 7.49862
122 // Par_2 = 14.7851 +/- 27.8012
124 // Parameterization Sy
126 // The parameterization of the error in the y direction along track uses
128 // #sigma_{y}^{||} = #sigma_{y}^{0} -a*exp(1/(x-b))
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
136 //==========================================================================
137 // Example how to retrive reference plots from the task
138 // void steerClErrParam(Int_t fig=0)
140 // gSystem->Load("libANALYSIS.so");
141 // gSystem->Load("libTRDqaRec.so");
143 // // initialize DB manager
144 // AliCDBManager *cdb = AliCDBManager::Instance();
145 // cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
147 // // initialize magnetic field.
148 // AliMagFCheb *field=new AliMagFCheb("Maps","Maps", 2, 1., 10., AliMagFCheb::k5kG);
149 // AliTracker::SetFieldMap(field, kTRUE);
151 // AliTRDclusterResolution *res = new AliTRDclusterResolution();
153 // res->Load("TRD.TaskClErrParam.root");
156 // //res->SetSaveAs();
157 // res->SetProcessCharge(kFALSE);
158 // res->SetProcessCenterPad(kFALSE);
159 // //res->SetProcessMean(kFALSE);
160 // res->SetProcessSigma(kFALSE);
161 // if(!res->PostProcess()) return;
163 // res->GetRefFigure(fig);
167 // Alexandru Bercuci <A.Bercuci@gsi.de> //
168 ////////////////////////////////////////////////////////////////////////////
170 #include "AliTRDclusterResolution.h"
171 #include "AliTRDtrackInfo/AliTRDclusterInfo.h"
172 #include "AliTRDgeometry.h"
173 #include "AliTRDcalibDB.h"
174 #include "AliTRDCommonParam.h"
175 #include "Cal/AliTRDCalROC.h"
176 #include "Cal/AliTRDCalDet.h"
179 #include "AliTracker.h"
180 #include "AliCDBManager.h"
183 #include "TObjArray.h"
186 #include "TGraphErrors.h"
190 #include "TLinearFitter.h"
196 ClassImp(AliTRDclusterResolution)
199 //_______________________________________________________
200 AliTRDclusterResolution::AliTRDclusterResolution(const char *name)
201 : AliTRDrecoTask(name, "Cluster Error Parametrization")
212 // ideal equidistant binning
213 //fAt = new TAxis(kNTB, -0.075, (kNTB-.5)*.15);
215 // equidistant binning for scaled for X0-MC (track ref)
216 fAt = new TAxis(kNTB, -0.075+.088, (kNTB-.5)*.15+.088);
218 // non-equidistant binning after vdrift correction
219 // const Double_t x[kNTB+1] = {
220 // 0.000000, 0.185084, 0.400490, 0.519701, 0.554653, 0.653150, 0.805063, 0.990261,
221 // 1.179610, 1.356406, 1.524094, 1.685499, 1.843083, 1.997338, 2.148077, 2.298274,
222 // 2.448656, 2.598639, 2.747809, 2.896596, 3.045380, 3.195000, 3.340303, 3.461688,
223 // 3.540094, 3.702677};
224 // fAt = new TAxis(kNTB, x);
226 fAd = new TAxis(kND, 0., .25);
228 // By default register all analysis
229 // The user can switch them off in his steering macro
231 SetProcessCenterPad();
236 //_______________________________________________________
237 AliTRDclusterResolution::~AliTRDclusterResolution()
239 if(fCanvas) delete fCanvas;
248 //_______________________________________________________
249 void AliTRDclusterResolution::ConnectInputData(Option_t *)
251 fInfo = dynamic_cast<TObjArray *>(GetInputData(0));
254 //_______________________________________________________
255 void AliTRDclusterResolution::CreateOutputObjects()
257 OpenFile(0, "RECREATE");
258 fContainer = Histos();
261 //_______________________________________________________
262 Bool_t AliTRDclusterResolution::GetRefFigure(Int_t ifig)
264 if(!fResults) return kFALSE;
267 TObjArray *arr = 0x0;
269 TGraphErrors *gm(0x0), *gs(0x0), *gp(0x0);
272 if(!(arr = (TObjArray*)fResults->At(kQRes))) break;
273 if(!(gm = (TGraphErrors*)arr->At(0))) break;
274 if(!(gs = (TGraphErrors*)arr->At(1))) break;
275 if(!(gp = (TGraphErrors*)arr->At(2))) break;
277 gs->GetHistogram()->GetYaxis()->SetRangeUser(-.01, .6);
278 gs->GetHistogram()->SetXTitle("Q [a.u.]");
279 gs->GetHistogram()->SetYTitle("#sigma_{y} / #mu_{y} [mm] / freq");
284 if(!(arr = (TObjArray*)fResults->At(kCenter))) break;
285 gPad->Divide(3, 1); l = gPad->GetListOfPrimitives();
286 for(Int_t ipad = 3; ipad--;){
287 if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
288 ((TVirtualPad*)l->At(ipad))->cd();
293 if(!(arr = (TObjArray*)fResults->At(kSigm))) break;
294 gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
295 for(Int_t ipad = 2; ipad--;){
296 if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
297 ((TVirtualPad*)l->At(ipad))->cd();
302 if(!(arr = (TObjArray*)fResults->At(kMean))) break;
303 gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
304 for(Int_t ipad = 2; ipad--;){
305 if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
306 ((TVirtualPad*)l->At(ipad))->cd();
313 AliWarning("No container/data found.");
317 //_______________________________________________________
318 TObjArray* AliTRDclusterResolution::Histos()
320 if(fContainer) return fContainer;
321 fContainer = new TObjArray(sizeof(AliTRDclusterResolution::EResultContainer));
322 //fContainer->SetOwner(kTRUE);
325 TObjArray *arr = 0x0;
327 fContainer->AddAt(h2 = new TH2I("h_q", "", 50, 2.2, 7.5, 100, -.5, .5), kQRes);
328 h2->SetXTitle("log(q) [a.u.]");
329 h2->SetYTitle("#Delta y[cm]");
330 h2->SetZTitle("entries");
332 fContainer->AddAt(arr = new TObjArray(AliTRDgeometry::kNlayer), kCenter);
333 arr->SetName("y(PadWidth)");
334 for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
335 arr->AddAt(h2 = new TH2I(Form("h_y%d", ily), Form("Ly[%d]", ily), 51, -.51, .51, 100, -.5, .5), ily);
336 h2->SetXTitle("y_{w} [w]");
337 h2->SetYTitle("#Delta y[cm]");
338 h2->SetZTitle("entries");
341 fContainer->AddAt(arr = new TObjArray(kN), kSigm);
342 arr->SetName("Resolution");
344 for(Int_t id=1; id<=fAd->GetNbins(); id++){
345 for(Int_t it=1; it<=fAt->GetNbins(); it++){
346 arr->AddAt(h2 = new TH2I(Form("h_d%02dt%02d", id, it), Form("d_{wire}(%3.1f-%3.1f)[mm] x_{drift}(%5.2f-%5.2f)[mm]", 10.*fAd->GetBinCenter(id)- 5.*fAd->GetBinWidth(id), 10.*fAd->GetBinCenter(id)+ 5.*fAd->GetBinWidth(id), 10.*fAt->GetBinCenter(it)- 5.*fAt->GetBinWidth(it), 10.*fAt->GetBinCenter(it)+ 5.*fAt->GetBinWidth(it)), 35, -.35, .35, 100, -.5, .5), ih++);
347 h2->SetXTitle("tg#phi");
348 h2->SetYTitle("#Delta y[cm]");
349 h2->SetZTitle("entries");
353 fContainer->AddAt(arr = new TObjArray(kN), kMean);
354 arr->SetName("Systematics");
356 for(Int_t id=1; id<=fAd->GetNbins(); id++){
357 for(Int_t it=1; it<=fAt->GetNbins(); it++){
358 arr->AddAt(h2 = new TH2I(Form("h_d%02dt%02d", id, it), Form("d_{wire}(%3.1f-%3.1f)[mm] x_{drift}(%5.2f-%5.2f)[mm]", 10.*fAd->GetBinCenter(id)- 5.*fAd->GetBinWidth(id), 10.*fAd->GetBinCenter(id)+ 5.*fAd->GetBinWidth(id), 10.*fAt->GetBinCenter(it)- 5.*fAt->GetBinWidth(it), 10.*fAt->GetBinCenter(it)+ 5.*fAt->GetBinWidth(it)), 35, -.35, .35, 100, -.5, .5), ih++);
359 h2->SetXTitle("tg#phi-h*tg(#theta)");
360 h2->SetYTitle("#Delta y[cm]");
361 h2->SetZTitle("entries");
368 //_______________________________________________________
369 void AliTRDclusterResolution::Exec(Option_t *)
371 if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the task.");
374 Float_t x, y, z, q, dy, dydx, dzdx, cov[3], covcl[3];
377 // define limits around ExB for which x contribution is negligible
378 const Float_t kDtgPhi = 3.5e-2; //(+- 2 deg)
380 TObjArray *arr0 = (TObjArray*)fContainer->At(kCenter);
381 TObjArray *arr1 = (TObjArray*)fContainer->At(kSigm);
382 TObjArray *arr2 = (TObjArray*)fContainer->At(kMean);
384 const AliTRDclusterInfo *cli = 0x0;
385 TIterator *iter=fInfo->MakeIterator();
386 while((cli=dynamic_cast<AliTRDclusterInfo*>((*iter)()))){
387 cli->GetCluster(det, x, y, z, q, t, covcl);
388 if(fDet>=0 && fDet!=det) continue;
390 dy = cli->GetResolution();
391 cli->GetGlobalPosition(y, z, dydx, dzdx, &cov[0]);
393 // resolution as a function of cluster charge
394 // only for phi equal exB
395 if(TMath::Abs(dydx-fExB) < kDtgPhi){
396 h2 = (TH2I*)fContainer->At(kQRes);
397 h2->Fill(TMath::Log(q), dy);
400 // do not use problematic clusters in resolution analysis
401 // TODO define limits as calibration aware (gain) !!
402 if(q<20. || q>250.) continue;
404 // resolution as a function of y displacement from pad center
405 // only for phi equal exB
406 if(TMath::Abs(dydx-fExB) < kDtgPhi){
407 h2 = (TH2I*)arr0->At(AliTRDgeometry::GetLayer(det));
408 h2->Fill(cli->GetYDisplacement(), dy);
411 Int_t it = fAt->FindBin(cli->GetDriftLength());
412 if(it==0 || it == fAt->GetNbins()+1){
413 AliWarning(Form("Drift length %f outside allowed range", cli->GetDriftLength()));
416 Int_t id = fAd->FindBin(cli->GetAnisochronity());
417 if(id==0 || id == fAd->GetNbins()+1){
418 AliWarning(Form("Distance to anode %f outside allowed range", cli->GetAnisochronity()));
421 // calculate index of cluster position in array
422 Int_t hid = (id-1)*kNTB+it-1;
424 // fill histo for resolution (sigma)
425 h2 = (TH2I*)arr1->At(hid);
428 // fill histo for systematic (mean)
429 h2 = (TH2I*)arr2->At(hid);
430 h2->Fill(dydx-cli->GetTilt()*dzdx, dy);
432 PostData(0, fContainer);
436 //_______________________________________________________
437 Bool_t AliTRDclusterResolution::PostProcess()
439 if(!fContainer) return kFALSE;
440 if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the post processing.");
442 TObjArray *arr = 0x0;
445 TGraphErrors *g = 0x0;
446 fResults = new TObjArray(sizeof(AliTRDclusterResolution::EResultContainer));
447 fResults->SetOwner();
448 fResults->AddAt(arr = new TObjArray(3), kQRes);
450 arr->AddAt(g = new TGraphErrors(), 0);
451 g->SetLineColor(kBlue); g->SetMarkerColor(kBlue);
452 g->SetMarkerStyle(7);
453 arr->AddAt(g = new TGraphErrors(), 1);
454 g->SetLineColor(kRed); g->SetMarkerColor(kRed);
455 g->SetMarkerStyle(23);
456 arr->AddAt(g = new TGraphErrors(), 2);
457 g->SetLineColor(kGreen); g->SetMarkerColor(kGreen);
458 g->SetMarkerStyle(7);
460 fResults->AddAt(arr = new TObjArray(3), kCenter);
463 if(!(h2 = (TH2F*)gROOT->FindObject("hYM"))){
464 h2 = new TH2F("hYM", "",
465 AliTRDgeometry::kNlayer, -.5, AliTRDgeometry::kNlayer-.5, 51, -.51, .51);
469 h2->SetYTitle("y [w]");
470 h2->SetZTitle("#mu_{x} [cm]");
471 arr->AddAt(h2 = (TH2F*)h2->Clone("hYS"), 1);
472 h2->SetZTitle("#sigma_{x} [cm]");
473 arr->AddAt(h2 = (TH2F*)h2->Clone("hYP"), 2);
474 h2->SetZTitle("entries");
476 fResults->AddAt(arr = new TObjArray(2), kSigm);
478 if(!(h2 = (TH2F*)gROOT->FindObject("hSX"))){
479 h2 = new TH2F("hSX", "",
480 fAd->GetNbins(), fAd->GetXmin(), fAd->GetXmax(),
481 fAt->GetNbins(), fAt->GetXmin(), fAt->GetXmax());
484 h2->SetXTitle("d [mm]");
485 h2->SetYTitle("x_{drift} [cm]");
486 h2->SetZTitle("#sigma_{x} [cm]");
487 arr->AddAt(h2 = (TH2F*)h2->Clone("hSY"), 1);
488 h2->SetZTitle("#sigma_{y} [cm]");
490 fResults->AddAt(arr = new TObjArray(2), kMean);
492 arr->AddAt(h2 = (TH2F*)h2->Clone("hDX"), 0);
493 h2->SetZTitle("dx [cm]");
494 arr->AddAt(h2 = (TH2F*)h2->Clone("hX0"), 1);
495 h2->SetZTitle("x0 [cm]");
498 TIterator *iter=fResults->MakeIterator();
499 while((o=((*iter)()))) o->Clear(); // maybe it is wrong but we should never reach this point
502 // precalculated value of tg^2(alpha_L)
503 Double_t exb2 = fExB*fExB;
504 // square of the mean value of sigma drift length.
505 // has to come from previous calibration
506 //Double_t sxd2 = 1.;// [mm^2]
508 printf("ExB[%e] ExB2[%e]\n", fExB, exb2);
510 // process resolution dependency on charge
511 if(HasProcessCharge()) ProcessCharge();
513 // process resolution dependency on y displacement
514 if(HasProcessCenterPad()) ProcessCenterPad();
516 // process resolution dependency on drift legth and drift cell width
517 if(HasProcessSigma()) ProcessSigma();
519 // process systematic shift on drift legth and drift cell width
520 if(HasProcessMean()) ProcessMean();
525 //_______________________________________________________
526 Bool_t AliTRDclusterResolution::SetExB(Int_t det, Int_t col, Int_t row)
529 AliCDBManager *cdb = AliCDBManager::Instance();
530 if(cdb->GetRun() < 0){
531 AliError("OCDB manager not properly initialized");
535 // check magnetic field
536 if(TMath::Abs(AliTracker::GetBz()) < 1.e-10){
537 AliWarning("B=0. Magnetic field may not be initialized. Continue if you know what you are doing !");
540 // set reference detector if any
541 if(det>=0 && det<AliTRDgeometry::kNdet) fDet = det;
544 AliTRDcalibDB *fCalibration = AliTRDcalibDB::Instance();
545 AliTRDCalROC *fCalVdriftROC = fCalibration->GetVdriftROC(det);
546 const AliTRDCalDet *fCalVdriftDet = fCalibration->GetVdriftDet();
548 fVdrift = fCalVdriftDet->GetValue(det) * fCalVdriftROC->GetValue(col, row);
549 fExB = AliTRDCommonParam::Instance()->GetOmegaTau(fVdrift);
554 //_______________________________________________________
555 void AliTRDclusterResolution::SetVisual()
558 fCanvas = new TCanvas("clResCanvas", "Cluster Resolution Visualization", 10, 10, 600, 600);
561 //_______________________________________________________
562 void AliTRDclusterResolution::ProcessCharge()
565 if((h2 = (TH2I*)fContainer->At(kQRes))) {
566 AliWarning("Missing dy=f(Q) histo");
569 TF1 f("f", "gaus", -.5, .5);
573 TObjArray *arr = (TObjArray*)fResults->At(kQRes);
574 TGraphErrors *gqm = (TGraphErrors*)arr->At(0);
575 TGraphErrors *gqs = (TGraphErrors*)arr->At(1);
576 TGraphErrors *gqp = (TGraphErrors*)arr->At(2);
577 Double_t q, n = 0., entries;
579 for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
580 Double_t q = TMath::Exp(ax->GetBinCenter(ix));
581 if(q<20. || q>250.) continue; // ?!
583 h1 = h2->ProjectionY("py", ix, ix);
584 entries = h1->GetEntries();
585 if(entries < 50) continue;
590 Int_t ip = gqm->GetN();
591 gqm->SetPoint(ip, q, 10.*f.GetParameter(1));
592 gqm->SetPointError(ip, 0., 10.*f.GetParError(1));
594 // correct sigma for ExB effect
595 gqs->SetPoint(ip, q, 1.e1*f.GetParameter(2)/**f.GetParameter(2)-exb2*sxd2*/);
596 gqs->SetPointError(ip, 0., 1.e1*f.GetParError(2)/**f.GetParameter(2)*/);
600 gqp->SetPoint(ip, q, entries);
601 gqp->SetPointError(ip, 0., 0./*TMath::Sqrt(entries)*/);
604 // normalize probability and get mean sy
605 Double_t sm = 0., sy;
606 for(Int_t ip=gqp->GetN(); ip--;){
607 gqp->GetPoint(ip, q, entries);
609 gqp->SetPoint(ip, q, entries);
610 gqs->GetPoint(ip, q, sy);
614 // error parametrization s(q) = <sy> + b(1/q-1/q0)
615 TF1 fq("fq", "[0] + [1]/x", 20., 250.);
617 printf("sy(Q) :: sm[%f] b[%f] 1/q0[%f]\n", sm, fq.GetParameter(1), (sm-fq.GetParameter(0))/fq.GetParameter(1));
620 //_______________________________________________________
621 void AliTRDclusterResolution::ProcessCenterPad()
623 TObjArray *arr = (TObjArray*)fContainer->At(kCenter);
625 AliWarning("Missing dy=f(y_d) container");
628 TF1 f("f", "gaus", -.5, .5);
630 TH1D *h1 = 0x0, *h11 = 0x0;
633 TObjArray *arrg = (TObjArray*)fResults->At(kCenter);
634 TH2F *hym = (TH2F*)arrg->At(0);
635 TH2F *hys = (TH2F*)arrg->At(1);
636 TH2F *hyp = (TH2F*)arrg->At(2);
637 for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
638 if(!(h2 = (TH2I*)arr->At(ily))) continue;
640 for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
641 Float_t yd = ax->GetBinCenter(ix);
642 h1 = h2->ProjectionY("py", ix, ix);
643 Int_t entries = (Int_t)h1->GetEntries();
644 if(entries < 50) continue;
649 hyp->Fill(ily, yd, entries);
650 hym->Fill(ily, yd, f.GetParameter(1));
651 //hym->SetPointError(ip, 0., f.GetParError(1));
652 hys->Fill(ily, yd, f.GetParameter(2));
653 //hys->SetPointError(ip, 0., f.GetParError(2));
657 // POSTPROCESS SPECTRA
658 // Found correction for systematic deviation
659 TF1 fprf("fprf", "[0]+[1]*sin([2]*x)", -.5, .5);
660 fprf.SetParameter(0, 0.);
661 fprf.SetParameter(1, 1.1E-2);
662 fprf.SetParameter(2, -TMath::PiOver2()/0.25);
663 printf(" const Float_t cy[AliTRDgeometry::kNlayer][3] = {\n");
664 for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
665 h1 = hym->ProjectionY("hym_pyy", ily+1, ily+1);
667 for(Int_t ib=h1->GetNbinsX(); ib--;) h1->SetBinError(ib, 0.002);
669 printf(" {%5.3e, %5.3e, %5.3e},\n", fprf.GetParameter(0), fprf.GetParameter(1), fprf.GetParameter(2));
671 if(!fCanvas) continue;
673 h1->SetMinimum(-0.02);h1->SetMaximum(0.02);h1->Draw("e1");
674 h11 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
675 h11->Scale(.8/h11->Integral());
676 h11->SetLineColor(kBlue); h11->Draw("csame");
677 fCanvas->Modified(); fCanvas->Update();
679 fCanvas->SaveAs(Form("Figures/ProcessCenterPad_M_Ly%d.gif", ily));
680 else gSystem->Sleep(500);
684 // Parameterization for sigma PRF
685 TF1 fgaus("fgaus", "gaus", -.5, .5);
686 printf(" const Float_t scy[AliTRDgeometry::kNlayer][4] = {\n");
687 for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
688 h1 = hys->ProjectionY("hys_pyy", ily+1, ily+1);
690 for(Int_t ib=h1->GetNbinsX(); ib--;) h1->SetBinError(ib, 0.002);
692 h1->Fit(&fgaus, "Q");
693 printf(" {%5.3e, %5.3e, %5.3e, ", fgaus.GetParameter(0), fgaus.GetParameter(1), fgaus.GetParameter(2));
695 // calculate mean sigma on the pad center distribution
697 h1 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
698 for(Int_t ix=1; ix<=h1->GetNbinsX(); ix++){
699 sy += fgaus.Eval(h1->GetBinCenter(ix))*h1->GetBinContent(ix);
701 sy /= h1->GetEntries();
702 printf("%5.3e},\n", sy);
704 if(!fCanvas) continue;
706 h1->SetMinimum(0.01);h1->SetMaximum(0.04);h1->Draw("e1");
707 h11 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
708 h11->Scale(1./h11->Integral());
709 h11->SetLineColor(kBlue); h11->Draw("csame");
710 fCanvas->Modified(); fCanvas->Update();
712 fCanvas->SaveAs(Form("Figures/ProcessCenterPad_S_Ly%d.gif", ily));
713 else gSystem->Sleep(500);
718 //_______________________________________________________
719 void AliTRDclusterResolution::ProcessSigma()
721 TObjArray *arr = (TObjArray*)fContainer->At(kSigm);
723 AliWarning("Missing dy=f(x_d, d_wire) container");
726 TLinearFitter gs(1,"pol1");
727 TF1 f("f", "gaus", -.5, .5);
731 // init visualization
732 TGraphErrors *ggs = 0x0;
735 ggs = new TGraphErrors();
739 Double_t d(0.), x(0.), sx, sy;
740 TObjArray *arrr = (TObjArray*)fResults->At(kSigm);
741 TH2F *hsx = (TH2F*)arrr->At(0);
742 TH2F *hsy = (TH2F*)arrr->At(1);
743 for(Int_t id=1; id<=fAd->GetNbins(); id++){
744 d = fAd->GetBinCenter(id); //[mm]
745 printf(" Doing d = %5.3f [mm]\n", d);
746 for(Int_t it=1; it<=fAt->GetNbins(); it++){
747 x = fAt->GetBinCenter(it); //[mm]
748 Int_t idx = (id-1)*kNTB+it-1;
750 // retrieve data histogram
751 if(!(h2 = (TH2I*)arr->At(idx))) {
752 AliWarning(Form("Missing histo at index idx[%3d] [id[%2d] it[%2d]] xd[%f] d[%f]\n", idx, id, it, x, d));
757 new(ggs) TGraphErrors();
758 ggs->SetMarkerStyle(7);
762 for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
763 Float_t dydx = ax->GetBinCenter(ix);
764 Double_t tgg = (dydx-fExB)/(1.+dydx*fExB);
765 Double_t tgg2 = tgg*tgg;
766 h1 = h2->ProjectionY("py", ix, ix);
767 if(h1->GetEntries() < 100) continue;
769 //printf("\tFit ix[%d] on %s entries [%d]\n", ix, h2->GetName(), (Int_t)h1->GetEntries());
771 Double_t s2 = f.GetParameter(2)*f.GetParameter(2);
772 Double_t s2e = 2.*f.GetParameter(2)*f.GetParError(2);
773 // Fill sy^2 = f(tg^2(phi-a_L))
774 gs.AddPoint(&tgg2, s2, s2e);
777 Int_t ip = ggs->GetN();
778 ggs->SetPoint(ip, tgg2, s2);
779 ggs->SetPointError(ip, 0., s2e);
780 //printf("%f, %f, %f, \n", tgg2, s2, s2e);
782 if(gs.Eval()) continue;
783 sx = gs.GetParameter(1);
784 sy = gs.GetParameter(0);
785 hsx->SetBinContent(id, it, TMath::Sqrt(sx));
786 //hsx->SetBinError(id, it, sig.GetParError(1));
788 // s^2_y = s0^2_y + tg^2(a_L) * s^2_x
789 hsy->SetBinContent(id, it, TMath::Sqrt(sy)/* - exb2*sig.GetParameter(1)*/);
790 //hsy->SetBinError(id, it, sig.GetParError(0)+exb2*exb2*sig.GetParError(1));
792 /*if(!fCanvas) */continue;
794 new(line) TLine(0, sy, .025, sy+.025*sx);
795 line->SetLineColor(kRed);line->SetLineWidth(2);
796 ggs->Draw("apl"); line->Draw("");
797 fCanvas->Modified(); fCanvas->Update();
798 if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessSigma_D%d_T%02d.gif", id, it));
799 else gSystem->Sleep(100);
801 printf(" xd=%4.1f[cm] sx=%5.3e[cm] sy=%5.3e[cm]\n", x, sx, sy);
804 // if(ggs) delete ggs; // TODO fix this memory leak !
805 // if(line) delete line;
807 // fit sy parallel to the drift
808 TF1 fsyD("fsy", "[0]+[1]*exp(1./(x+[2]))", 0.5, 3.5);
809 fsyD.SetParameter(0, .025);
810 fsyD.SetParameter(1, -8.e-4);
811 fsyD.SetParameter(2, -.25);
812 h1 = hsy->ProjectionY("hsy_pyy"); h1->Scale(1./hsy->GetNbinsX());
813 for(Int_t ib=1; ib<=h1->GetNbinsX(); ib++) h1->SetBinError(ib, 0.002);
814 h1->Fit(&fsyD, "Q", "", .5, 3.5);
815 printf(" const Float_t sy0 = %5.3e;\n", fsyD.GetParameter(0));
816 printf(" const Float_t sya = %5.3e;\n", fsyD.GetParameter(1));
817 printf(" const Float_t syb = %5.3e;\n", fsyD.GetParameter(2));
820 h1->Draw("e1"); fsyD.Draw("same");
821 fCanvas->Modified(); fCanvas->Update();
822 if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SY||.gif");
823 else gSystem->Sleep(1000);
827 // fit sx parallel to the drift
828 TF1 fsxD("fsx", "gaus(0)+expo(3)", 0., 3.5);
829 h1 = hsx->ProjectionY("hsx_pyy");
830 h1->Scale(1./hsx->GetNbinsX());
831 for(Int_t ib=1; ib<=h1->GetNbinsX(); ib++) if(h1->GetBinContent(ib)>0.) h1->SetBinError(ib, 0.02);
832 h1->Fit(&f, "QN", "", 0., 1.3);
833 fsxD.SetParameter(0, 0.);
834 fsxD.SetParameter(1, f.GetParameter(1));
835 fsxD.SetParameter(2, f.GetParameter(2));
836 TF1 expo("fexpo", "expo", 0., 3.5);
837 h1->Fit(&expo, "QN");
838 fsxD.SetParameter(3, expo.GetParameter(0));
839 fsxD.SetParameter(4, expo.GetParameter(1));
841 printf(" const Float_t sxgc = %5.3e;\n", fsxD.GetParameter(0));
842 printf(" const Float_t sxgm = %5.3e;\n", fsxD.GetParameter(1));
843 printf(" const Float_t sxgs = %5.3e;\n", fsxD.GetParameter(2));
844 printf(" const Float_t sxe0 = %5.3e;\n", fsxD.GetParameter(3));
845 printf(" const Float_t sxe1 = %5.3e;\n", fsxD.GetParameter(4));
848 h1->Draw("e1"); fsxD.Draw("same");
849 fCanvas->Modified(); fCanvas->Update();
850 if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SX||.gif");
851 else gSystem->Sleep(1000);
854 // fit sx perpendicular to the drift
856 TAxis *ay = hsx->GetYaxis();
857 TH1D *hx = hsx->ProjectionX("hsx_px", 1,1);
858 TH1D *hAn = (TH1D*)hx->Clone("hAn"); hAn->Reset();
859 for(Int_t ib = 11; ib<24; ib++, nb++){
860 hx = hsx->ProjectionX("hsx_px", ib, ib);
861 for(Int_t id=1; id<=hx->GetNbinsX(); id++)
862 hx->SetBinContent(id, hx->GetBinContent(id)-fsxD.Eval(ay->GetBinCenter(ib)));
865 TF1 fsxA("fsxA", "pol2", 0., 0.25);
867 for(Int_t ib=1; ib<=hAn->GetNbinsX(); ib++) hAn->SetBinError(ib, 0.002);
868 hAn->Fit(&fsxA, "Q");
869 printf(" const Float_t sxd0 = %5.3e;\n", fsxA.GetParameter(0));
870 printf(" const Float_t sxd1 = %5.3e;\n", fsxA.GetParameter(1));
871 printf(" const Float_t sxd2 = %5.3e;\n", fsxA.GetParameter(2));
874 hAn->Draw("e1"); fsxA.Draw("same");
875 fCanvas->Modified(); fCanvas->Update();
876 if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SXL.gif");
877 else gSystem->Sleep(100);
881 //_______________________________________________________
882 void AliTRDclusterResolution::ProcessMean()
884 TObjArray *arr = (TObjArray*)fContainer->At(kMean);
886 AliWarning("Missing dy=f(x_d, d_wire) container");
889 TGraphErrors *gm = new TGraphErrors();
890 TF1 f("f", "gaus", -.5, .5);
891 TF1 line("l", Form("%6.4f*[0]+[1]*x", fExB), -.15, .15);
892 Double_t d(0.), x(0.), dx, x0;
897 TObjArray *arrr = (TObjArray*)fResults->At(kMean);
898 TH2F *hdx = (TH2F*)arrr->At(0);
899 TH2F *hx0 = (TH2F*)arrr->At(1);
900 for(Int_t id=1; id<=fAd->GetNbins(); id++){
901 d = fAd->GetBinCenter(id); //[mm]
902 printf(" Doing d = %5.3f [mm]\n", d);
903 for(Int_t it=1; it<=fAt->GetNbins(); it++){
904 x = fAt->GetBinCenter(it); //[mm]
905 Int_t idx = (id-1)*kNTB+it-1;
906 if(!(h2 = (TH2I*)arr->At(idx))) {
907 AliWarning(Form("Missing histo at index idx[%3d] [id[%2d] it[%2d]] xd[%f] d[%f]\n", idx, id, it, x, d));
911 new(gm) TGraphErrors();
912 gm->SetMarkerStyle(7);
914 for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
915 Double_t dydx = ax->GetBinCenter(ix);
916 h1 = h2->ProjectionY("py", ix, ix);
917 if(h1->GetEntries() < 50) continue;
921 // Fill <Dy> = f(dydx - h*dzdx)
922 Int_t ip = gm->GetN();
923 gm->SetPoint(ip, dydx, f.GetParameter(1));
924 gm->SetPointError(ip, 0., f.GetParError(1));
927 if(gm->GetN()<4) continue;
929 gm->Fit(&line, "QN");
930 dx = line.GetParameter(1); // = (fVdrift + dVd)*(fT + tCorr);
931 x0 = line.GetParameter(0); // = tg(a_L)*dx
932 hdx->SetBinContent(id, it, dx);
933 hx0->SetBinContent(id, it, x0);
935 if(!fCanvas) continue;
937 gm->Draw("apl"); line.Draw("same");
938 fCanvas->Modified(); fCanvas->Update();
939 if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessMean_D%d_T%02d.gif", id, it));
940 else gSystem->Sleep(100);
941 printf(" xd=%4.1f[cm] dx=%5.3e[cm] x0=%5.3e[cm]\n", x, dx, x0);
945 TH1D *h=hdx->ProjectionY("hdx_py"); h->Scale(1./hdx->GetNbinsX());
946 printf(" const Float_t cx[] = {\n ");
947 for(Int_t ib=1; ib<=h->GetNbinsX(); ib++){
948 printf("%6.4e, ", h->GetBinContent(ib));
949 if(ib%6==0) printf("\n ");
953 // delete gm; TODO memory leak ?