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