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