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