]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDclusterResolution.cxx
copy TRD performance train to PWG1
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDclusterResolution.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercialf purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 /* $Id: 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"
173 #include "AliTRDcluster.h"
174 #include "AliTRDcalibDB.h"
175 #include "AliTRDCommonParam.h"
176 #include "Cal/AliTRDCalROC.h"
177 #include "Cal/AliTRDCalDet.h"
178
179 #include "AliLog.h"
180 #include "AliTracker.h"
181 #include "AliCDBManager.h"
182
183 #include "TROOT.h"
184 #include "TObjArray.h"
185 #include "TAxis.h"
186 #include "TF1.h"
187 #include "TLegend.h"
188 #include "TGraphErrors.h"
189 #include "TLine.h"
190 #include "TH2I.h"
191 #include "TH3S.h"
192 #include "TTree.h"
193 #include "TMath.h"
194 #include "TLinearFitter.h"
195
196 #include "TCanvas.h"
197 #include "TSystem.h"
198
199 ClassImp(AliTRDclusterResolution)
200
201 const Float_t AliTRDclusterResolution::fgkTimeBinLength = 1./ AliTRDCommonParam::Instance()->GetSamplingFrequency();
202 //_______________________________________________________
203 AliTRDclusterResolution::AliTRDclusterResolution(const char *name, const char *title)
204   : AliTRDrecoTask(name, title)
205   ,fCanvas(0x0)
206   ,fInfo(0x0)
207   ,fResults(0x0)
208   ,fAt(0x0)
209   ,fStatus(0)
210   ,fDet(-1)
211   ,fExB(0.)
212   ,fVdrift(0.)
213   ,fLy(0)
214   ,fX(0.)
215   ,fY(0.)
216   ,fZ(0.)
217 {
218 // Constructor
219
220   memset(fR, 0, 4*sizeof(Float_t));
221   memset(fP, 0, 4*sizeof(Float_t));
222   // time drift axis
223   fAt = new TAxis(kNTB, 0., kNTB*fgkTimeBinLength);
224
225   // By default register all analysis
226   // The user can switch them off in his steering macro
227   SetProcess(kQRes);
228   SetProcess(kCenter);
229   SetProcess(kMean);
230   SetProcess(kSigm);
231 }
232
233 //_______________________________________________________
234 AliTRDclusterResolution::~AliTRDclusterResolution()
235 {
236 // Destructor
237
238   if(fCanvas) delete fCanvas;
239   if(fAt) delete fAt;
240   if(fResults){
241     fResults->Delete();
242     delete fResults;
243   }
244 }
245
246 //_______________________________________________________
247 void AliTRDclusterResolution::ConnectInputData(Option_t *)
248 {
249   fInfo = dynamic_cast<TObjArray *>(GetInputData(0));
250 }
251
252 //_______________________________________________________
253 void AliTRDclusterResolution::CreateOutputObjects()
254 {
255   OpenFile(0, "RECREATE");
256   fContainer = Histos();
257 }
258
259 //_______________________________________________________
260 Bool_t AliTRDclusterResolution::GetRefFigure(Int_t ifig)
261 {
262 // Steering function to retrieve performance plots
263
264   if(!fResults) return kFALSE;
265   TLegend *leg = 0x0;
266   TList *l = 0x0;
267   TObjArray *arr = 0x0;
268   TTree *t = 0x0;
269   TH2 *h2 = 0x0;TH1 *h1 = 0x0;
270   TGraphErrors *gm(0x0), *gs(0x0), *gp(0x0);
271   switch(ifig){
272   case kQRes:
273     if(!(arr = (TObjArray*)fResults->At(kQRes))) break;
274     if(!(gm = (TGraphErrors*)arr->At(0))) break;
275     if(!(gs = (TGraphErrors*)arr->At(1))) break;
276     if(!(gp = (TGraphErrors*)arr->At(2))) break;
277     gs->Draw("apl");
278     gs->GetHistogram()->GetYaxis()->SetRangeUser(-50., 700.);
279     gs->GetHistogram()->SetXTitle("Q [a.u.]");
280     gs->GetHistogram()->SetYTitle("#sigma_{y} / #mu_{y} [#mum] / freq");
281     gm->Draw("pl");
282     gp->Draw("pl");
283     return kTRUE;
284   case kCenter:
285     if(!(arr = (TObjArray*)fResults->At(kCenter))) break;
286     gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
287     ((TVirtualPad*)l->At(0))->cd();
288     ((TTree*)arr->At(0))->Draw("y:x>>h(23, 0.1, 2.4, 51, -.51, .51)",
289             "m[0]*(ly==0&&abs(m[0])<1.e-1)", "colz");
290     ((TVirtualPad*)l->At(1))->cd();
291     leg= new TLegend(.7, .7, .9, .95);
292     leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetFillStyle(0);
293     leg->SetHeader("TRD Plane"); 
294     for(Int_t il = 1; il<=AliTRDgeometry::kNlayer; il++){
295       if(!(gm = (TGraphErrors*)arr->At(il))) return kFALSE;
296       gm->Draw(il>1?"pc":"apc"); leg->AddEntry(gm, Form("%d", il-1), "pl");
297       if(il>1) continue;
298       gm->GetHistogram()->SetXTitle("t_{drift} [#mus]");
299       gm->GetHistogram()->SetYTitle("#sigma_{y}(x|cen=0) [#mum]");
300       gm->GetHistogram()->GetYaxis()->SetRangeUser(150., 500.);
301     }
302     leg->Draw();
303     return kTRUE;
304   case kSigm:
305     if(!(t = (TTree*)fResults->At(kSigm))) break;
306     t->Draw("z:t>>h2x(23, 0.1, 2.4, 25, 0., 2.5)","sx*(1)", "lego2fb");
307     h2 = (TH2F*)gROOT->FindObject("h2x");
308     printf("  const Double_t sx[24][25]={\n");
309     for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
310       printf("    {");
311       for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
312         printf("%6.4f ", h2->GetBinContent(ix, iy));
313       }
314       printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
315     }
316     printf("  };\n");
317     gPad->Divide(2, 1, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
318     ((TVirtualPad*)l->At(0))->cd();
319     h1 = h2->ProjectionX("hsx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
320     h1->SetYTitle("<#sigma_{x}> [#mum]");
321     h1->SetXTitle("t_{drift} [#mus]");
322     h1->GetXaxis()->SetRange(2, kNTB-1); h1->Draw("pc");
323
324     t->Draw("z:t>>h2y(23, 0.1, 2.4, 25, 0., 2.5)","sy*(1)", "lego2fb");
325     h2 = (TH2F*)gROOT->FindObject("h2y");
326     printf("  const Double_t sy[24][25]={\n");
327     for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
328       printf("    {");
329       for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
330         printf("%6.4f ", h2->GetBinContent(ix, iy));
331       }
332       printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
333     }
334     printf("  };\n");
335     ((TVirtualPad*)l->At(1))->cd();
336     h1 = h2->ProjectionX("hsy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
337     h1->SetYTitle("<#sigma_{y}> [#mum]");
338     h1->SetXTitle("t_{drift} [#mus]");
339     h1->GetXaxis()->SetRange(2, kNTB-1); h1->Draw("pc");
340     return kTRUE;
341   case kMean:
342     if(!(t = (TTree*)fResults->At(kMean))) break;
343     t->Draw("z:t>>h2x(23, 0.1, 2.4, 25, 0., 2.5)","dx*(1)", "goff");
344     h2 = (TH2F*)gROOT->FindObject("h2x");
345     printf("  const Double_t dx[24][25]={\n");
346     for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
347       printf("    {");
348       for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
349         printf("%6.4f ", h2->GetBinContent(ix, iy));
350       }
351       printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
352     }
353     printf("  };\n");
354     gPad->Divide(2, 1, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
355     ((TVirtualPad*)l->At(0))->cd();
356     h1 = h2->ProjectionX("hdx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
357     h1->SetYTitle("<dx> [#mum]");
358     h1->SetXTitle("t_{drift} [#mus]");
359     h1->GetXaxis()->SetRange(2, kNTB-1); h1->Draw("pc");
360
361     t->Draw("z:t>>h2y(23, 0.1, 2.4, 25, 0., 2.5)","dy*(1)", "goff");
362     h2 = (TH2F*)gROOT->FindObject("h2y");
363     printf("  const Double_t dy[24][25]={\n");
364     for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
365       printf("    {");
366       for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
367         printf("%6.4f ", h2->GetBinContent(ix, iy));
368       }
369       printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
370     }
371     printf("  };\n");
372     ((TVirtualPad*)l->At(1))->cd();
373     h1 = h2->ProjectionX("hdy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
374     h1->SetYTitle("<dy> [#mum]");
375     h1->SetXTitle("t_{drift} [#mus]");
376     h1->GetXaxis()->SetRange(2, kNTB-1); h1->Draw("pc");
377
378     return kTRUE;
379   default:
380     break;
381   }
382   AliWarning("No container/data found.");
383   return kFALSE;
384 }
385
386 //_______________________________________________________
387 TObjArray* AliTRDclusterResolution::Histos()
388 {
389 // Retrieve histograms array if already build or build it
390
391   if(fContainer) return fContainer;
392   fContainer = new TObjArray(kNtasks);
393   //fContainer->SetOwner(kTRUE);
394
395   TH3S *h3 = 0x0;
396   TObjArray *arr = 0x0;
397
398   fContainer->AddAt(arr = new TObjArray(2*AliTRDgeometry::kNlayer), kCenter);
399   arr->SetName("Center");
400   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
401     // add resolution plot for each layer
402     if(!(h3=(TH3S*)gROOT->FindObject(Form("hCenResLy%d", il)))){ 
403       h3 = new TH3S(
404         Form("hCenResLy%d", il), 
405         Form(" ly [%d]", il), 
406         kNTB, fAt->GetBinLowEdge(1), fAt->GetBinUpEdge(kNTB),   // x
407         51, -.51, .51, // y 
408         60, -.3, .3); // dy
409       h3->SetXTitle("x [#mus]");
410       h3->SetYTitle("y [pw]");
411       h3->SetZTitle("#Delta y[cm]");
412     } h3->Reset();
413     arr->AddAt(h3, il);
414     // add Pull plot for each layer
415     if(!(h3=(TH3S*)gROOT->FindObject(Form("hCenPullLy%d", il)))){ 
416       h3 = new TH3S(
417         Form("hCenPullLy%d", il), 
418         Form(" ly [%d]", il), 
419         kNTB, fAt->GetBinLowEdge(1), fAt->GetBinUpEdge(kNTB),   // x
420         51, -.51, .51, // y 
421         60, -4., 4.); // dy
422       h3->SetXTitle("x [#mus]");
423       h3->SetYTitle("y [pw]");
424       h3->SetZTitle("#Delta y/#sigma_{y}");
425     } h3->Reset();
426     arr->AddAt(h3, AliTRDgeometry::kNlayer+il);
427   }
428
429   if(!(h3 = (TH3S*)gROOT->FindObject("Charge"))){
430     h3 = new TH3S("Charge", "dy=f(q)", 50, 2.2, 7.5, 60, -.3, .3, 60, -4., 4.);
431     h3->SetXTitle("log(q) [a.u.]");
432     h3->SetYTitle("#Delta y[cm]");
433     h3->SetZTitle("#Delta y/#sigma_{y}");
434   }
435   fContainer->AddAt(h3, kQRes);
436
437   fContainer->AddAt(arr = new TObjArray(kNTB), kSigm);
438   arr->SetName("Resolution");
439   for(Int_t ix=0; ix<kNTB; ix++){
440     if(!(h3=(TH3S*)gROOT->FindObject(Form("hr_x%02d", ix)))){
441       h3 = new TH3S(
442         Form("hr_x%02d", ix), 
443         Form(" t_{drift}(%3.1f-%3.1f)[#mus]", fAt->GetBinLowEdge(ix+1), fAt->GetBinUpEdge(ix+1)), 
444         kND, 0., 2.5,   // z 
445         35, -.35, .35, // tgp 
446         60, -.3, .3); // dy
447       h3->SetXTitle("z [mm]");
448       h3->SetYTitle("tg#phi");
449       h3->SetZTitle("#Delta y[cm]");
450     }
451     arr->AddAt(h3, ix);
452   }
453
454   fContainer->AddAt(arr = new TObjArray(kNTB), kMean);
455   arr->SetName("Systematics");
456   for(Int_t ix=0; ix<kNTB; ix++){
457     if(!(h3=(TH3S*)gROOT->FindObject(Form("hs_x%02d", ix)))){
458       h3 = new TH3S(
459         Form("hs_x%02d", ix), 
460         Form(" t_{drift}(%3.1f-%3.1f)[#mus]", fAt->GetBinLowEdge(ix+1), fAt->GetBinUpEdge(ix+1)), 
461         kND, 0., 2.5,   // z 
462         35, -.35, .35, // tgp-h tgt 
463         60, -.3, .3); // dy
464       h3->SetXTitle("z [mm]");
465       h3->SetYTitle("tg(#phi) - h*tg(#theta)");
466       h3->SetZTitle("#Delta y[cm]");
467     }
468     arr->AddAt(h3, ix);
469   }
470
471   return fContainer;
472 }
473
474 //_______________________________________________________
475 void AliTRDclusterResolution::Exec(Option_t *)
476 {
477 // Fill container histograms
478
479   if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the task.");
480
481   Int_t det, t;
482   Float_t x, y, z, q, dy, dydx, dzdx, cov[3], covcl[3];
483   TH3S *h3 = 0x0;
484
485   // define limits around ExB for which x contribution is negligible
486   const Float_t kDtgPhi = 3.5e-2; //(+- 2 deg)
487
488   TObjArray *arr0 = (TObjArray*)fContainer->At(kCenter);
489   TObjArray *arr1 = (TObjArray*)fContainer->At(kSigm);
490   TObjArray *arr2 = (TObjArray*)fContainer->At(kMean);
491
492   const AliTRDclusterInfo *cli = 0x0;
493   TIterator *iter=fInfo->MakeIterator();
494   while((cli=dynamic_cast<AliTRDclusterInfo*>((*iter)()))){
495     cli->GetCluster(det, x, y, z, q, t, covcl);
496     if(fDet>=0 && fDet!=det) continue;
497     
498     dy = cli->GetResolution();
499     cli->GetGlobalPosition(y, z, dydx, dzdx, &cov[0]);
500
501     // resolution as a function of cluster charge
502     // only for phi equal exB 
503     if(TMath::Abs(dydx-fExB) < kDtgPhi){
504       h3 = (TH3S*)fContainer->At(kQRes);
505       h3->Fill(TMath::Log(q), dy, dy/TMath::Sqrt(covcl[0]));
506
507       printf("q=%f Log(q)=%f dy=%f pull=%f\n",q, TMath::Log(q), dy, dy/TMath::Sqrt(covcl[0]));
508     }
509
510     // do not use problematic clusters in resolution analysis
511     // TODO define limits as calibration aware (gain) !!
512     if(q<20. || q>250.) continue;
513
514     x = (t+.5)*fgkTimeBinLength; // conservative approach !!
515
516     // resolution as a function of y displacement from pad center
517     // only for phi equal exB
518     if(TMath::Abs(dydx-fExB) < kDtgPhi/* &&
519        TMath::Abs(x-0.675)<0.225*/){
520       Int_t ly(AliTRDgeometry::GetLayer(det));
521       h3 = (TH3S*)arr0->At(ly);
522       h3->Fill(x, cli->GetYDisplacement(), dy);
523       h3 = (TH3S*)arr0->At(AliTRDgeometry::kNlayer+ly);
524       h3->Fill(x, cli->GetYDisplacement(), dy/TMath::Sqrt(covcl[0]));
525     }
526
527     Int_t ix = fAt->FindBin(x);
528     if(ix==0 || ix == fAt->GetNbins()+1){
529       AliWarning(Form("Drift time %3.1f outside allowed range", x));
530       continue;
531     }
532
533     // fill histo for resolution (sigma)
534     ((TH3S*)arr1->At(ix-1))->Fill(10.*cli->GetAnisochronity(), dydx, dy);
535
536     // fill histo for systematic (mean)
537     ((TH3S*)arr2->At(ix-1))->Fill(10.*cli->GetAnisochronity(), dydx-cli->GetTilt()*dzdx, dy);  
538   }
539   PostData(0, fContainer);
540 }
541
542
543 //_______________________________________________________
544 Bool_t AliTRDclusterResolution::PostProcess()
545 {
546   if(!fContainer) return kFALSE;
547   if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the post processing.");
548   
549   TObjArray *arr = 0x0;
550   TTree *t=0x0;
551   if(!fResults){
552     TGraphErrors *g = 0x0;
553     fResults = new TObjArray(kNtasks);
554     fResults->SetOwner();
555     fResults->AddAt(arr = new TObjArray(3), kQRes);
556     arr->SetOwner();
557     arr->AddAt(g = new TGraphErrors(), 0);
558     g->SetLineColor(kBlue); g->SetMarkerColor(kBlue);
559     g->SetMarkerStyle(7); 
560     arr->AddAt(g = new TGraphErrors(), 1);
561     g->SetLineColor(kRed); g->SetMarkerColor(kRed);
562     g->SetMarkerStyle(23); 
563     arr->AddAt(g = new TGraphErrors(), 2);
564     g->SetLineColor(kGreen); g->SetMarkerColor(kGreen);
565     g->SetMarkerStyle(7); 
566
567     // pad center dependence
568     fResults->AddAt(arr = new TObjArray(AliTRDgeometry::kNlayer+1), kCenter);
569     arr->SetOwner();
570     arr->AddAt(
571     t = new TTree("cent", "dy=f(y,x,ly)"), 0);
572     t->Branch("ly", &fLy, "ly/B");
573     t->Branch("x", &fX, "x/F");
574     t->Branch("y", &fY, "y/F");
575     t->Branch("m", &fR[0], "m[2]/F");
576     t->Branch("s", &fR[2], "s[2]/F");
577     t->Branch("pm", &fP[0], "pm[2]/F");
578     t->Branch("ps", &fP[2], "ps[2]/F");
579     for(Int_t il=1; il<=AliTRDgeometry::kNlayer; il++){
580       arr->AddAt(g = new TGraphErrors(), il);
581       g->SetLineColor(il); g->SetLineStyle(il);
582       g->SetMarkerColor(il);g->SetMarkerStyle(4); 
583     }
584
585
586     fResults->AddAt(t = new TTree("sigm", "dy=f(dw,x,dydx)"), kSigm);
587     t->Branch("t", &fX, "t/F");
588     t->Branch("z", &fZ, "z/F");
589     t->Branch("sx", &fR[0], "sx[2]/F");
590     t->Branch("sy", &fR[2], "sy[2]/F");
591
592
593     fResults->AddAt(t = new TTree("mean", "dy=f(dw,x,dydx - h dzdx)"), kMean);
594     t->Branch("t", &fX, "t/F");
595     t->Branch("z", &fZ, "z/F");
596     t->Branch("dx", &fR[0], "dx[2]/F");
597     t->Branch("dy", &fR[2], "dy[2]/F");
598   } else {
599     TObject *o = 0x0;
600     TIterator *iter=fResults->MakeIterator();
601     while((o=((*iter)()))) o->Clear(); // maybe it is wrong but we should never reach this point
602   }
603   
604   // precalculated value of tg^2(alpha_L)
605   Double_t exb2 = fExB*fExB;
606   // square of the mean value of sigma drift length.
607   // has to come from previous calibration 
608   //Double_t sxd2 = 1.;// [mm^2]
609
610   printf("ExB[%e] ExB2[%e]\n", fExB, exb2);
611
612   // process resolution dependency on charge
613   if(HasProcess(kQRes)) ProcessCharge();
614   
615   // process resolution dependency on y displacement
616   if(HasProcess(kCenter)) ProcessCenterPad();
617
618   // process resolution dependency on drift legth and drift cell width
619   if(HasProcess(kSigm)) ProcessSigma();
620
621   // process systematic shift on drift legth and drift cell width
622   if(HasProcess(kMean)) ProcessMean();
623
624   return kTRUE;
625 }
626
627 //_______________________________________________________
628 Bool_t AliTRDclusterResolution::SetExB(Int_t det, Int_t col, Int_t row)
629 {
630   // check OCDB
631   AliCDBManager *cdb = AliCDBManager::Instance();
632   if(cdb->GetRun() < 0){
633     AliError("OCDB manager not properly initialized");
634     return kFALSE;
635   }
636
637   // check magnetic field
638   if(TMath::Abs(AliTracker::GetBz()) < 1.e-10){
639     AliWarning("B=0. Magnetic field may not be initialized. Continue if you know what you are doing !");
640   }
641
642   // set reference detector if any
643   if(det>=0 && det<AliTRDgeometry::kNdet) fDet = det;
644   else det = 0;
645
646   AliTRDcalibDB *fCalibration  = AliTRDcalibDB::Instance();
647   AliTRDCalROC  *fCalVdriftROC = fCalibration->GetVdriftROC(det);
648   const AliTRDCalDet  *fCalVdriftDet = fCalibration->GetVdriftDet();
649
650   fVdrift = fCalVdriftDet->GetValue(det) * fCalVdriftROC->GetValue(col, row);
651   fExB   = AliTRDCommonParam::Instance()->GetOmegaTau(fVdrift);
652   SetBit(kExB);
653   return kTRUE;
654 }
655
656 //_______________________________________________________
657 void AliTRDclusterResolution::SetVisual()
658 {
659   if(fCanvas) return;
660   fCanvas = new TCanvas("clResCanvas", "Cluster Resolution Visualization", 10, 10, 600, 600);
661 }
662
663 //_______________________________________________________
664 void AliTRDclusterResolution::ProcessCharge()
665 {
666 // Resolution as a function of cluster charge.
667 //
668 // As described in the function ProcessCenterPad() the error parameterization for clusters for phi = a_L can be 
669 // written as:
670 // BEGIN_LATEX
671 // #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
672 // END_LATEX
673 // with the contribution in case of B=0 given by:
674 // BEGIN_LATEX
675 // #sigma_{y}|_{B=0} = #sigma_{diff}*Gauss(0, s_{ly}) + #delta_{#sigma}(q)
676 // END_LATEX
677 // which further can be simplified to:
678 // BEGIN_LATEX
679 // <#sigma_{y}|_{B=0}>(q) = <#sigma_{y}> + #delta_{#sigma}(q)
680 // <#sigma_{y}> = #int{f(q)#sigma_{y}dq}
681 // END_LATEX
682 // The results for s_y and f(q) are displayed below:
683 //Begin_Html
684 //<img src="TRD/clusterQerror.gif">
685 //End_Html
686 // The function has to extended to accomodate gain calibration scalling and errors.
687 //
688 // Author
689 // Alexandru Bercuci <A.Bercuci@gsi.de>
690
691   TH2I *h2 = 0x0;
692   if(!(h2 = (TH2I*)fContainer->At(kQRes))) {
693     AliWarning("Missing dy=f(Q) histo");
694     return;
695   }
696   TF1 f("f", "gaus", -.5, .5);
697   TAxis *ax = 0x0;
698   TH1D *h1 = 0x0;
699
700   // compute mean error on x
701   Double_t s2x = 0.; 
702   for(Int_t ix=5; ix<kNTB; ix++){
703     // retrieve error on the drift length
704     s2x += AliTRDcluster::GetSX(ix);
705   }
706   s2x /= (kNTB-5); s2x *= s2x;
707   Double_t exb2 = fExB*fExB;
708
709   TObjArray *arr = (TObjArray*)fResults->At(kQRes);
710   TGraphErrors *gqm = (TGraphErrors*)arr->At(0);
711   TGraphErrors *gqs = (TGraphErrors*)arr->At(1);
712   TGraphErrors *gqp = (TGraphErrors*)arr->At(2);
713   Double_t q, n = 0., entries;
714   ax = h2->GetXaxis();
715   for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
716     q = TMath::Exp(ax->GetBinCenter(ix));
717     if(q<20. || q>250.) continue; // ?!
718
719     h1 = h2->ProjectionY("py", ix, ix);
720     entries = h1->GetEntries();
721     if(entries < 50) continue;
722     Adjust(&f, h1);
723     h1->Fit(&f, "Q");
724
725     // Fill sy^2 = f(q)
726     Int_t ip = gqm->GetN();
727     gqm->SetPoint(ip, q, 1.e4*f.GetParameter(1));
728     gqm->SetPointError(ip, 0., 1.e4*f.GetParError(1));
729
730     // correct sigma for ExB effect
731     gqs->SetPoint(ip, q, 1.e4*(f.GetParameter(2)*f.GetParameter(2)-exb2*s2x));
732     gqs->SetPointError(ip, 0., 1.e4*f.GetParError(2)*f.GetParameter(2));
733
734     // save probability
735     n += entries;
736     gqp->SetPoint(ip, q, entries);
737     gqp->SetPointError(ip, 0., 0./*TMath::Sqrt(entries)*/);
738   } 
739
740   // normalize probability and get mean sy
741   Double_t sm = 0., sy;
742   for(Int_t ip=gqp->GetN(); ip--;){
743     gqp->GetPoint(ip, q, entries);
744     entries/=n;
745     gqp->SetPoint(ip, q, 1.e3*entries);
746     gqs->GetPoint(ip, q, sy);
747     sm += entries*sy;
748   }
749
750   // error parametrization s(q) = <sy> + b(1/q-1/q0)
751   TF1 fq("fq", "[0] + [1]/x", 20., 250.);
752   gqs->Fit(&fq/*, "W"*/);
753   printf("sm=%f [0]=%f [1]=%f\n", 1.e-4*sm, fq.GetParameter(0), fq.GetParameter(1));
754   printf("  const Float_t sq0inv = %f; // [1/q0]\n", (sm-fq.GetParameter(0))/fq.GetParameter(1));
755   printf("  const Float_t sqb    = %f; // [cm]\n", 1.e-4*fq.GetParameter(1));
756 }
757
758 //_______________________________________________________
759 void AliTRDclusterResolution::ProcessCenterPad()
760 {
761 // Resolution as a function of y displacement from pad center and drift length.
762 //
763 // Since the error parameterization of cluster r-phi position can be written as (see AliTRDcluster::SetSigmaY2()):
764 // BEGIN_LATEX
765 // #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
766 // END_LATEX
767 // one can see that for phi = a_L one gets the following expression:
768 // BEGIN_LATEX
769 // #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
770 // END_LATEX
771 // where we have explicitely marked the remaining term in case of absence of magnetic field. Thus one can use the 
772 // previous equation to estimate s_y for B=0 and than by comparing in magnetic field conditions one can get the s_x.
773 // This is a simplified method to determine the error parameterization for s_x and s_y as compared to the one 
774 // implemented in ProcessSigma(). For more details on cluster error parameterization please see also 
775 // AliTRDcluster::SetSigmaY2()
776 // 
777 // The representation of dy=f(y_cen, x_drift| layer) can be also used to estimate the systematic shift in the r-phi 
778 // coordinate resulting from imperfection in the cluster shape parameterization. From the expresion of the shift derived 
779 // in ProcessMean() with phi=exb one gets: 
780 // BEGIN_LATEX
781 // <#Delta y>= <#delta x> * (tg(#alpha_{L})-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
782 // <#Delta y>(y_{cen})= -h*<#delta x>(x_{drift}, q_{cl}) * dz/dx + #delta y(y_{cen}, ...)
783 // END_LATEX
784 // where all dependences are made explicit. This last expression can be used in two ways:
785 //   - by average on the dz/dx we can determine directly dy (the method implemented here) 
786 //   - by plotting as a function of dzdx one can determine both dx and dy components in an independent method.
787 //Begin_Html
788 //<img src="TRD/clusterYcorr.gif">
789 //End_Html
790 // Author
791 // Alexandru Bercuci <A.Bercuci@gsi.de>
792
793   TObjArray *arr = (TObjArray*)fContainer->At(kCenter);
794   if(!arr) {
795     AliWarning("Missing dy=f(y | x, ly) container");
796     return;
797   }
798   Double_t exb2 = fExB*fExB;
799   Float_t s[AliTRDgeometry::kNlayer];
800   TF1 f("f", "gaus", -.5, .5);
801   TF1 fp("fp", "gaus", -3.5, 3.5);
802
803   TH1D *h1 = 0x0; TH2F *h2 = 0x0; TH3S *h3r=0x0, *h3p=0x0;
804   TObjArray *arrRes = (TObjArray*)fResults->At(kCenter);
805   TTree *t = (TTree*)arrRes->At(0);
806   TGraphErrors *gs = 0x0;
807   TAxis *ax = 0x0;
808
809   printf("  const Float_t lSy[6][24] = {\n      {");
810   const Int_t nl = AliTRDgeometry::kNlayer;
811   for(Int_t il=0; il<nl; il++){
812     if(!(h3r = (TH3S*)arr->At(il))) continue;
813     if(!(h3p = (TH3S*)arr->At(nl+il))) continue;
814     gs = (TGraphErrors*)arrRes->At(il+1);
815     fLy = il;
816 //    printf("Ly[%d]\n", il);
817     for(Int_t ix=1; ix<=h3r->GetXaxis()->GetNbins(); ix++){
818       ax = h3r->GetXaxis(); ax->SetRange(ix, ix);
819       ax = h3p->GetXaxis(); ax->SetRange(ix, ix);
820       fX  = ax->GetBinCenter(ix);
821 //      printf("  x[%2d]=%4.2f\n", ix, fX);
822       for(Int_t iy=1; iy<=h3r->GetYaxis()->GetNbins(); iy++){
823         ax = h3r->GetYaxis(); ax->SetRange(iy, iy);
824         ax = h3p->GetYaxis(); ax->SetRange(iy, iy);
825         fY  = ax->GetBinCenter(iy);
826 //        printf("    y[%2d]=%5.2f\n", iy, fY);
827         // finish navigation in the HnSparse
828
829         h1 = (TH1D*)h3r->Project3D("z");
830         Int_t entries = (Int_t)h1->Integral();
831         if(entries < 50) continue;
832         //Adjust(&f, h1);
833         h1->Fit(&f, "QN");
834     
835         // Fill sy,my=f(y_w,x,ly)
836         fR[0] = f.GetParameter(1); fR[1] = f.GetParError(1);
837         fR[2] = f.GetParameter(2); fR[3] = f.GetParError(2);
838
839         h1 = (TH1D*)h3p->Project3D("z");
840         h1->Fit(&fp, "QN");
841         fP[0] = fp.GetParameter(1); fP[1] = fp.GetParError(1);
842         fP[2] = fp.GetParameter(2); fP[3] = fp.GetParError(2);
843
844         //printf("ly[%d] x[%3.1f] y[%+5.2f] m[%5.3f] s[%5.3f] \n", fLy, fX, fY, fR[0], fR[2]);
845         t->Fill();
846
847
848       }
849     }
850     t->Draw("y:x>>h(24, 0., 2.4, 51, -.51, .51)",
851             Form("s[0]*(ly==%d&&abs(m[0])<1.e-1)", fLy),
852             "goff");
853     h2=(TH2F*)gROOT->FindObject("h");
854     f.FixParameter(1, 0.);
855     Int_t n = h2->GetXaxis()->GetNbins(), nn(0); s[il]=0.;
856     printf("    {");
857     for(Int_t ix=1; ix<=n; ix++){
858       ax = h2->GetXaxis();
859       fX  = ax->GetBinCenter(ix);
860       h1 = h2->ProjectionY("hCenPy", ix, ix);
861       //if((Int_t)h1->Integral() < 1.e-10) continue; 
862
863       // Apply lorentz angle correction
864       // retrieve error on the drift length
865       Double_t s2x = AliTRDcluster::GetSX(ix-1); s2x *= s2x;
866       Int_t nnn = 0;
867       for(Int_t iy=1; iy<=h1->GetNbinsX(); iy++){
868         Double_t s2 = h1->GetBinContent(iy); s2*= s2;
869         // sigma square corrected for Lorentz angle
870         // s2 = s2_y(y_w,x)+exb2*s2_x
871         Double_t sy = TMath::Sqrt(TMath::Max(s2 - exb2*s2x, Double_t(0.)));
872         if(sy<1.e-20) continue;
873         h1->SetBinContent(iy, sy); nnn++;
874         printf("s[%6.2f] sx[%6.2f] sy[%6.2f]\n",
875         1.e4*TMath::Sqrt(s2), 1.e4*TMath::Abs(fExB*AliTRDcluster::GetSX(ix-1)), 
876         1.e4*h1->GetBinContent(iy));
877       }
878       // do fit only if enough data
879       Double_t sPRF = 0.;
880       if(nnn>5){
881         h1->Fit(&f, "QN");
882         sPRF = f.GetParameter(2);
883         nn++;
884       }
885       s[il]+=sPRF;
886       printf("%6.4f,%s", sPRF, ix%6?" ":"\n     ");
887       Int_t jx = gs->GetN();
888       gs->SetPoint(jx, fX, 1.e4*sPRF);
889       gs->SetPointError(jx, 0., 0./*f.GetParError(0)*/);
890     }
891     printf("\b},\n");
892     s[il]/=nn;
893
894     f.ReleaseParameter(2);
895
896
897     if(!fCanvas) continue;
898     h2->Draw("lego2fb");
899     fCanvas->Modified(); fCanvas->Update();
900     if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessCenter_ly[%d].gif", fLy));
901     else gSystem->Sleep(100);
902   }
903   printf("  };\n");
904   printf("  const Float_t lPRF[] = {"
905     "%5.3f, %5.3f, %5.3f, %5.3f, %5.3f, %5.3f};\n",
906     s[0], s[1], s[2], s[3], s[4], s[5]);
907 }
908
909 //_______________________________________________________
910 void AliTRDclusterResolution::ProcessSigma()
911 {
912 // As the r-phi coordinate is the only one which is measured by the TRD detector we have to rely on it to
913 // estimate both the radial (x) and r-phi (y) errors. This method is based on the following assumptions. 
914 // The measured error in the y direction is the sum of the intrinsic contribution of the r-phi measurement
915 // with the contribution of the radial measurement - because x is not a parameter of Alice track model (Kalman).
916 // BEGIN_LATEX
917 // #sigma^{2}|_{y} = #sigma^{2}_{y*} + #sigma^{2}_{x*}   
918 // END_LATEX
919 // In the general case 
920 // BEGIN_LATEX
921 // #sigma^{2}_{y*} = #sigma^{2}_{y} + tg^{2}(#alpha_{L})#sigma^{2}_{x_{drift}}   
922 // #sigma^{2}_{x*} = tg^{2}(#phi - #alpha_{L})*(#sigma^{2}_{x_{drift}} + #sigma^{2}_{x_{0}} + tg^{2}(#alpha_{L})*x^{2}/12)
923 // END_LATEX
924 // where we have explicitely show the lorentz angle correction on y and the projection of radial component on the y
925 // direction through the track angle in the bending plane (phi). Also we have shown that the radial component in the
926 // last equation has twp terms, the drift and the misalignment (x_0). For ideal geometry or known misalignment one 
927 // can solve the equation
928 // BEGIN_LATEX
929 // #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}]
930 // END_LATEX
931 // by fitting a straight line:
932 // BEGIN_LATEX
933 // #sigma^{2}|_{y} = a(x_{cl}, z_{cl}) * tg^{2}(#phi - #alpha_{L}) + b(x_{cl}, z_{cl})
934 // END_LATEX
935 // the error parameterization will be given by:
936 // BEGIN_LATEX
937 // #sigma_{x} (x_{cl}, z_{cl}) = #sqrt{a(x_{cl}, z_{cl}) - tg^{2}(#alpha_{L})*x^{2}/12}
938 // #sigma_{y} (x_{cl}, z_{cl}) = #sqrt{b(x_{cl}, z_{cl}) - #sigma^{2}_{x} (x_{cl}, z_{cl}) * tg^{2}(#alpha_{L})}
939 // END_LATEX
940 // Below there is an example of such dependency. 
941 //Begin_Html
942 //<img src="TRD/clusterSigmaMethod.gif">
943 //End_Html
944 //
945 // The error parameterization obtained by this method are implemented in the functions AliTRDcluster::GetSX() and
946 // AliTRDcluster::GetSYdrift(). For an independent method to determine s_y as a function of drift length check the 
947 // function ProcessCenterPad(). One has to keep in mind that while this method return the mean s_y over the distance
948 // to pad center distribution the other method returns the *STANDARD* value at center=0 (maximum). To recover the 
949 // standard value one has to solve the obvious equation:
950 // BEGIN_LATEX
951 // #sigma_{y}^{STANDARD} = #frac{<#sigma_{y}>}{#int{s exp(s^{2}/#sigma) ds}}
952 // END_LATEX
953 // with "<s_y>" being the value calculated here and "sigma" the width of the s_y distribution calculated in 
954 // ProcessCenterPad().
955 //  
956 // Author
957 // Alexandru Bercuci <A.Bercuci@gsi.de>
958
959   TObjArray *arr = (TObjArray*)fContainer->At(kSigm);
960   if(!arr){
961     AliWarning("Missing dy=f(x_d, d_w) container");
962     return;
963   }
964
965   // init visualization
966   TGraphErrors *ggs = 0x0;
967   TGraph *line = 0x0;
968   if(fCanvas){
969     ggs = new TGraphErrors();
970     line = new TGraph();
971     line->SetLineColor(kRed);line->SetLineWidth(2);
972   }
973
974   // init logistic support
975   TF1 f("f", "gaus", -.5, .5);
976   TLinearFitter gs(1,"pol1");
977   TH1 *hFrame=0x0;
978   TH1D *h1 = 0x0; TH3S *h3=0x0;
979   TAxis *ax = 0x0;
980   Double_t exb2 = fExB*fExB, x;
981   AliTRDcluster c;
982   TTree *t = (TTree*)fResults->At(kSigm);
983   for(Int_t ix=0; ix<kNTB; ix++){
984     if(!(h3=(TH3S*)arr->At(ix))) continue;
985     c.SetPadTime(ix);
986     x = c.GetXloc(0., 1.5);
987     fX= fAt->GetBinCenter(ix+1);
988     for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
989       ax = h3->GetXaxis();
990       ax->SetRange(iz, iz);
991       fZ = ax->GetBinCenter(iz);
992
993       // reset visualization
994       if(fCanvas){ 
995         new(ggs) TGraphErrors();
996         ggs->SetMarkerStyle(7);
997       }
998       gs.ClearPoints();
999
1000       for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
1001         ax = h3->GetYaxis();
1002         ax->SetRange(ip, ip); 
1003         Double_t tgl = ax->GetBinCenter(ip);
1004         // finish navigation in the HnSparse
1005
1006         //if(TMath::Abs(dydx)>0.18) continue;
1007         Double_t tgg = (tgl-fExB)/(1.+tgl*fExB);
1008         Double_t tgg2 = tgg*tgg;
1009
1010         h1 = (TH1D*)h3->Project3D("z");
1011         Int_t entries = (Int_t)h1->Integral();
1012         if(entries < 50) continue;
1013         //Adjust(&f, h1);
1014         h1->Fit(&f, "QN");
1015
1016         Double_t s2  = f.GetParameter(2)*f.GetParameter(2);
1017         Double_t s2e = 2.*f.GetParameter(2)*f.GetParError(2);
1018         // Fill sy^2 = f(tg^2(phi-a_L))
1019         gs.AddPoint(&tgg2, s2, s2e);
1020
1021         if(!ggs) continue;
1022         Int_t jp = ggs->GetN();
1023         ggs->SetPoint(jp, tgg2, s2);
1024         ggs->SetPointError(jp, 0., s2e);
1025       }
1026       // TODO here a more robust fit method has to be provided
1027       // for which lower boundaries on the parameters have to 
1028       // be imposed. Unfortunately the Minuit fit does not work 
1029       // for the TGraph in the case of B not 0.
1030       if(gs.Eval()) continue;
1031
1032       fR[0] = gs.GetParameter(1) - x*x*exb2/12.;
1033       printf("s2x+x2=%f ang=%f s2x=%f\n", gs.GetParameter(1), x*x*exb2/12., fR[0]);
1034       fR[0] = TMath::Max(fR[0], Float_t(4.e-4)); 
1035
1036       // s^2_y  = s0^2_y + tg^2(a_L) * s^2_x
1037       // s0^2_y = f(D_L)*x + s_PRF^2 
1038       fR[2]= gs.GetParameter(0)-exb2*fR[0];
1039       printf("s2y+s2x=%f s2y=%f\n", fR[0], fR[2]);
1040       fR[2] = TMath::Max(fR[2], Float_t(2.5e-5)); 
1041       fR[0] = TMath::Sqrt(fR[0]);
1042       fR[1] = .5*gs.GetParError(1)/fR[0];
1043       fR[2] = TMath::Sqrt(fR[2]);
1044       fR[3] = gs.GetParError(0)+exb2*exb2*gs.GetParError(1);
1045       t->Fill();
1046       printf("    xd=%4.2f[cm] sx=%6.1f[um] sy=%5.1f[um]\n", x, 1.e4*fR[0], 1.e4*fR[2]);
1047
1048       if(!fCanvas) continue;
1049       fCanvas->cd(); fCanvas->SetLogx(); //fCanvas->SetLogy();
1050       if(!hFrame){ 
1051         fCanvas->SetMargin(0.15, 0.01, 0.1, 0.01);
1052         hFrame=new TH1I("hFrame", "", 100, 0., .3);
1053         hFrame->SetMinimum(0.);hFrame->SetMaximum(.005);
1054         hFrame->SetXTitle("tg^{2}(#phi-#alpha_{L})");
1055         hFrame->SetYTitle("#sigma^{2}y[cm^{2}]");
1056         hFrame->GetYaxis()->SetTitleOffset(2.);
1057         hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
1058         hFrame->Draw();
1059       } else hFrame->Reset();
1060       Double_t xx = 0., dxx=.2/50;
1061       for(Int_t ip=0;ip<50;ip++){ 
1062         line->SetPoint(ip, xx,  gs.GetParameter(0)+xx*gs.GetParameter(1)); 
1063         xx+=dxx;
1064       }
1065       ggs->Draw("pl"); line->Draw("l");
1066       fCanvas->Modified(); fCanvas->Update();
1067       if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessSigma_z[%5.3f]_x[%5.3f].gif", fZ, fX));
1068       else gSystem->Sleep(100);
1069     }
1070   }
1071   return;
1072 }
1073
1074 //_______________________________________________________
1075 void AliTRDclusterResolution::ProcessMean()
1076 {
1077 // By this method the cluster shift in r-phi and radial directions can be estimated by comparing with the MC.
1078 // The resolution of the cluster corrected for pad tilt with respect to MC in the r-phi (measuring) plane can be 
1079 // expressed by:
1080 // BEGIN_LATEX
1081 // #Delta y=w - y_{MC}(x_{cl})
1082 // w = y_{cl}^{'} + h*(z_{MC}(x_{cl})-z_{cl})
1083 // y_{MC}(x_{cl}) = y_{0} - dy/dx*x_{cl}
1084 // z_{MC}(x_{cl}) = z_{0} - dz/dx*x_{cl}
1085 // y_{cl}^{'} = y_{cl}-x_{cl}*tg(#alpha_{L})
1086 // END_LATEX
1087 // where x_cl is the drift length attached to a cluster, y_cl is the r-phi coordinate of the cluster measured by
1088 // charge sharing on adjacent pads and y_0 and z_0 are MC reference points (as example the track references at 
1089 // entrance/exit of a chamber). If we suppose that both r-phi (y) and radial (x) coordinate of the clusters are 
1090 // affected by errors we can write
1091 // BEGIN_LATEX
1092 // x_{cl} = x_{cl}^{*} + #delta x 
1093 // y_{cl} = y_{cl}^{*} + #delta y 
1094 // END_LATEX 
1095 // where the starred components are the corrected values. Thus by definition the following quantity
1096 // BEGIN_LATEX
1097 // #Delta y^{*}= w^{*} - y_{MC}(x_{cl}^{*})
1098 // END_LATEX
1099 // has 0 average over all dependency. Using this decomposition we can write:
1100 // BEGIN_LATEX
1101 // <#Delta y>=<#Delta y^{*}> + <#delta x * (dy/dx-h*dz/dx) + #delta y - #delta x * tg(#alpha_{L})>
1102 // END_LATEX
1103 // which can be transformed to the following linear dependence:
1104 // BEGIN_LATEX
1105 // <#Delta y>= <#delta x> * (dy/dx-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
1106 // END_LATEX
1107 // if expressed as function of dy/dx-h*dz/dx. Furtheremore this expression can be plotted for various clusters
1108 // i.e. we can explicitely introduce the diffusion (x_cl) and drift cell - anisochronity (z_cl) dependences. From 
1109 // plotting this dependence and linear fitting it with:
1110 // BEGIN_LATEX
1111 // <#Delta y>= a(x_{cl}, z_{cl}) * (dy/dx-h*dz/dx) + b(x_{cl}, z_{cl})
1112 // END_LATEX
1113 // the systematic shifts will be given by:
1114 // BEGIN_LATEX
1115 // #delta x (x_{cl}, z_{cl}) = a(x_{cl}, z_{cl})
1116 // #delta y (x_{cl}, z_{cl}) = b(x_{cl}, z_{cl}) + a(x_{cl}, z_{cl}) * tg(#alpha_{L})
1117 // END_LATEX
1118 // Below there is an example of such dependency. 
1119 //Begin_Html
1120 //<img src="TRD/clusterShiftMethod.gif">
1121 //End_Html
1122 //
1123 // The occurance of the radial shift is due to the following conditions 
1124 //   - the approximation of a constant drift velocity over the drift length (larger drift velocities close to 
1125 //     cathode wire plane)
1126 //   - the superposition of charge tails in the amplification region (first clusters appear to be located at the 
1127 //     anode wire)
1128 //   - the superposition of charge tails in the drift region (shift towards anode wire)
1129 //   - diffusion effects which convolute with the TRF thus enlarging it
1130 //   - approximate knowledge of the TRF (approximate measuring in test beam conditions) 
1131 // 
1132 // The occurance of the r-phi shift is due to the following conditions 
1133 //   - approximate model for cluster shape (LUT)
1134 //   - rounding-up problems
1135 //
1136 // The numerical results for ideal simulations for the radial and r-phi shifts are displayed below and used 
1137 // for the cluster reconstruction (see the functions AliTRDcluster::GetXcorr() and AliTRDcluster::GetYcorr()). 
1138 //Begin_Html
1139 //<img src="TRD/clusterShiftX.gif">
1140 //<img src="TRD/clusterShiftY.gif">
1141 //End_Html
1142 // More details can be found in the presentation given during the TRD
1143 // software meeting at the end of 2008 and beginning of year 2009, published on indico.cern.ch.
1144 // 
1145 // Author 
1146 // Alexandru Bercuci <A.Bercuci@gsi.de>
1147
1148
1149  
1150   TObjArray *arr = (TObjArray*)fContainer->At(kMean);
1151   if(!arr){
1152     AliWarning("Missing dy=f(x_d, d_w) container");
1153     return;
1154   }
1155
1156   // init logistic support
1157   TF1 f("f", "gaus", -.5, .5);
1158   TF1 line("l", "[0]+[1]*x", -.15, .15);
1159   TGraphErrors *gm = new TGraphErrors();
1160   TH1 *hFrame=0x0;
1161   TH1D *h1 = 0x0; TH3S *h3 =0x0;
1162   TAxis *ax = 0x0;
1163   Double_t x;
1164
1165   AliTRDcluster c;
1166   TTree *t = (TTree*)fResults->At(kMean);
1167   for(Int_t ix=0; ix<kNTB; ix++){
1168     if(!(h3=(TH3S*)arr->At(ix))) continue;
1169     c.SetPadTime(ix);
1170     x = c.GetXloc(0., 1.5);
1171     fX= fAt->GetBinCenter(ix+1);
1172     for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
1173       ax = h3->GetXaxis();
1174       ax->SetRange(iz, iz);
1175       fZ = ax->GetBinCenter(iz);
1176
1177       // reset fitter
1178       new(gm) TGraphErrors();
1179       gm->SetMarkerStyle(7);
1180
1181       for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
1182         ax = h3->GetYaxis();
1183         ax->SetRange(ip, ip); 
1184         Double_t tgl = ax->GetBinCenter(ip);
1185         // finish navigation in the HnSparse
1186
1187         h1 = (TH1D*)h3->Project3D("z");
1188         Int_t entries = (Int_t)h1->Integral();
1189         if(entries < 50) continue;
1190         //Adjust(&f, h1);
1191         h1->Fit(&f, "QN");
1192
1193         // Fill <Dy> = f(dydx - h*dzdx)
1194         Int_t jp = gm->GetN();
1195         gm->SetPoint(jp, tgl, f.GetParameter(1));
1196         gm->SetPointError(jp, 0., f.GetParError(1));
1197       }
1198       if(gm->GetN()<4) continue;
1199
1200       gm->Fit(&line, "QN");
1201       fR[0] = line.GetParameter(1); // dx
1202       fR[1] = line.GetParError(1);
1203       fR[2] = line.GetParameter(0) + fExB*fR[0]; // xs = dy - tg(a_L)*dx
1204       t->Fill();
1205       printf("    xd=%4.2f[cm] dx=%6.2f[um] dy=%6.2f[um]\n", x, 1.e4*fR[0], 1.e4*fR[2]);
1206
1207       if(!fCanvas) continue;
1208       fCanvas->cd();
1209       if(!hFrame){ 
1210         fCanvas->SetMargin(0.1, 0.02, 0.1, 0.01);
1211         hFrame=new TH1I("hFrame", "", 100, -.3, .3);
1212         hFrame->SetMinimum(-.1);hFrame->SetMaximum(.1);
1213         hFrame->SetXTitle("tg#phi-htg#theta");
1214         hFrame->SetYTitle("#Delta y[cm]");
1215         hFrame->GetYaxis()->SetTitleOffset(1.5);
1216         hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
1217         hFrame->Draw();
1218       } else hFrame->Reset();
1219       gm->Draw("pl"); line.Draw("same");
1220       fCanvas->Modified(); fCanvas->Update();
1221       if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessMean_Z[%5.3f]_X[%5.3f].gif", fZ, fX));
1222       else gSystem->Sleep(100);
1223     }
1224   }
1225 }