]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDclusterResolution.cxx
557d854ad999a5bf85e756d75e3f997175a2aa25
[u/mrichter/AliRoot.git] / TRD / qaRec / 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");
146 //   cdb->SetRun(0);
147 //   // initialize magnetic field.
148 //   AliMagFMaps *field = new AliMagFMaps("Maps","Maps", 2, 1., 10., AliMagFMaps::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 //   if(!res->PostProcess()) return;
156 //   new TCanvas;
157 //   res->GetRefFigure(fig);
158 // }
159 //
160 //  Authors:                                                              //
161 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
162 ////////////////////////////////////////////////////////////////////////////
163
164 #include "AliTRDclusterResolution.h"
165 #include "AliTRDtrackInfo/AliTRDclusterInfo.h"
166 #include "AliTRDgeometry.h"
167 #include "AliTRDcalibDB.h"
168 #include "Cal/AliTRDCalROC.h"
169 #include "Cal/AliTRDCalDet.h"
170
171 #include "AliLog.h"
172 #include "AliTracker.h"
173 #include "AliCDBManager.h"
174
175 #include "TROOT.h"
176 #include "TObjArray.h"
177 #include "TAxis.h"
178 #include "TF1.h"
179 #include "TGraphErrors.h"
180 #include "TLine.h"
181 #include "TH2I.h"
182 #include "TMath.h"
183 #include "TLinearFitter.h"
184
185 #include "TCanvas.h"
186 #include "TSystem.h"
187
188
189 ClassImp(AliTRDclusterResolution)
190
191
192 //_______________________________________________________
193 AliTRDclusterResolution::AliTRDclusterResolution()
194   : AliTRDrecoTask("ClErrParam", "Cluster Error Parametrization")
195   ,fCanvas(0x0)
196   ,fInfo(0x0)
197   ,fResults(0x0)
198   ,fAt(0x0)
199   ,fAd(0x0)
200   ,fStatus(0)
201   ,fDet(-1)
202   ,fExB(0.)
203   ,fVdrift(0.)
204 {
205   // ideal equidistant binning
206   //fAt = new TAxis(kNTB, -0.075, (kNTB-.5)*.15);
207
208   // equidistant binning for scaled for X0-MC (track ref)
209   fAt = new TAxis(kNTB, -0.075+.088, (kNTB-.5)*.15+.088);
210   
211   // non-equidistant binning after vdrift correction
212 //   const Double_t x[kNTB+1] = {
213 // 0.000000, 0.185084, 0.400490, 0.519701, 0.554653, 0.653150, 0.805063, 0.990261,
214 // 1.179610, 1.356406, 1.524094, 1.685499, 1.843083, 1.997338, 2.148077, 2.298274,
215 // 2.448656, 2.598639, 2.747809, 2.896596, 3.045380, 3.195000, 3.340303, 3.461688,
216 // 3.540094, 3.702677};
217 //   fAt = new TAxis(kNTB, x);
218
219   fAd = new TAxis(kND, 0., .25);
220
221   // By default register all analysis
222   // The user can switch them off in his steering macro
223   SetProcessCharge();
224   SetProcessCenterPad();
225   SetProcessMean();
226   SetProcessSigma();
227 }
228
229 //_______________________________________________________
230 AliTRDclusterResolution::~AliTRDclusterResolution()
231 {
232   if(fCanvas) delete fCanvas;
233   if(fAd) delete fAd;
234   if(fAt) delete fAt;
235   if(fResults){
236     fResults->Delete();
237     delete fResults;
238   }
239 }
240
241 //_______________________________________________________
242 void AliTRDclusterResolution::ConnectInputData(Option_t *)
243 {
244   fInfo = dynamic_cast<TObjArray *>(GetInputData(0));
245 }
246
247 //_______________________________________________________
248 void AliTRDclusterResolution::CreateOutputObjects()
249 {
250   OpenFile(0, "RECREATE");
251   fContainer = Histos();
252 }
253
254 //_______________________________________________________
255 Bool_t AliTRDclusterResolution::GetRefFigure(Int_t ifig)
256 {
257   if(!fResults) return kFALSE;
258   
259   TList *l = 0x0;
260   TObjArray *arr = 0x0;
261   TH2 *h2 = 0x0;
262   TGraphErrors *gm(0x0), *gs(0x0), *gp(0x0);
263   switch(ifig){
264   case kQRes:
265     if(!(arr = (TObjArray*)fResults->At(kQRes))) break;
266     if(!(gm = (TGraphErrors*)arr->At(0))) break;
267     if(!(gs = (TGraphErrors*)arr->At(1))) break;
268     if(!(gp = (TGraphErrors*)arr->At(2))) break;
269     gs->Draw("apl");
270     gs->GetHistogram()->GetYaxis()->SetRangeUser(-.01, .6);
271     gs->GetHistogram()->SetXTitle("Q [a.u.]");
272     gs->GetHistogram()->SetYTitle("#sigma_{y} / #mu_{y} [mm] / freq");
273     gm->Draw("pl");
274     gp->Draw("pl");
275     return kTRUE;
276   case kCenter:
277     if(!(arr = (TObjArray*)fResults->At(kCenter))) break;
278     gPad->Divide(3, 1); l = gPad->GetListOfPrimitives();
279     for(Int_t ipad = 3; ipad--;){
280       if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
281       ((TVirtualPad*)l->At(ipad))->cd();
282       h2->Draw("lego2fb");
283     }
284     return kTRUE;
285   case kSigm:
286     if(!(arr = (TObjArray*)fResults->At(kSigm))) break;
287     gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
288     for(Int_t ipad = 2; ipad--;){
289       if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
290       ((TVirtualPad*)l->At(ipad))->cd();
291       h2->Draw("lego2fb");
292     }
293     return kTRUE;
294   case kMean:
295     if(!(arr = (TObjArray*)fResults->At(kMean))) break;
296     gPad->Divide(2, 1);  l = gPad->GetListOfPrimitives();
297     for(Int_t ipad = 2; ipad--;){
298       if(!(h2 = (TH2F*)arr->At(ipad))) return kFALSE;
299       ((TVirtualPad*)l->At(ipad))->cd();
300       h2->Draw("lego2fb");
301     }
302     return kTRUE;
303   default:
304     break;
305   }
306   AliWarning("No container/data found.");
307   return kFALSE;
308 }
309
310 //_______________________________________________________
311 TObjArray* AliTRDclusterResolution::Histos()
312 {
313   if(fContainer) return fContainer;
314   fContainer = new TObjArray(sizeof(AliTRDclusterResolution::EResultContainers));
315   //fContainer->SetOwner(kTRUE);
316
317   TH2I *h2 = 0x0;
318   TObjArray *arr = 0x0;
319
320   fContainer->AddAt(h2 = new TH2I("h_q", "", 50, 2.2, 7.5, 100, -.5, .5), kQRes);
321   h2->SetXTitle("log(q) [a.u.]");
322   h2->SetYTitle("#Delta y[cm]");
323   h2->SetZTitle("entries");
324
325   fContainer->AddAt(arr = new TObjArray(AliTRDgeometry::kNlayer), kCenter);
326   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
327     arr->AddAt(h2 = new TH2I(Form("h_y%d", ily), Form("Ly[%d]", ily), 51, -.51, .51, 100, -.5, .5), ily);
328     h2->SetXTitle("y_{w} [w]");
329     h2->SetYTitle("#Delta y[cm]");
330     h2->SetZTitle("entries");
331   }
332
333   fContainer->AddAt(arr = new TObjArray(kN), kSigm);
334   Int_t ih = 0;
335   for(Int_t id=1; id<=fAd->GetNbins(); id++){
336     for(Int_t it=1; it<=fAt->GetNbins(); it++){
337       arr->AddAt(h2 = new TH2I(Form("h_d%02dt%02d", id, it), Form("d_{wire}(%3.1f-%3.1f)[mm] x_{drift}(%5.2f-%5.2f)[mm]", 10.*fAd->GetBinCenter(id)- 5.*fAd->GetBinWidth(id), 10.*fAd->GetBinCenter(id)+ 5.*fAd->GetBinWidth(id), 10.*fAt->GetBinCenter(it)- 5.*fAt->GetBinWidth(it), 10.*fAt->GetBinCenter(it)+ 5.*fAt->GetBinWidth(it)), 35, -.35, .35, 100, -.5, .5), ih++);
338       h2->SetXTitle("tg#phi");
339       h2->SetYTitle("#Delta y[cm]");
340       h2->SetZTitle("entries");
341     }
342   }
343
344   fContainer->AddAt(arr = new TObjArray(kN), kMean);
345   ih = 0;
346   for(Int_t id=1; id<=fAd->GetNbins(); id++){
347     for(Int_t it=1; it<=fAt->GetNbins(); it++){
348       arr->AddAt(h2 = new TH2I(Form("h_d%02dt%02d", id, it), Form("d_{wire}(%3.1f-%3.1f)[mm] x_{drift}(%5.2f-%5.2f)[mm]", 10.*fAd->GetBinCenter(id)- 5.*fAd->GetBinWidth(id), 10.*fAd->GetBinCenter(id)+ 5.*fAd->GetBinWidth(id), 10.*fAt->GetBinCenter(it)- 5.*fAt->GetBinWidth(it), 10.*fAt->GetBinCenter(it)+ 5.*fAt->GetBinWidth(it)), 35, -.35, .35, 100, -.5, .5), ih++);
349       h2->SetXTitle("tg#phi-h*tg(#theta)");
350       h2->SetYTitle("#Delta y[cm]");
351       h2->SetZTitle("entries");
352     }
353   }
354
355   return fContainer;
356 }
357
358 //_______________________________________________________
359 void AliTRDclusterResolution::Exec(Option_t *)
360 {
361   if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the task.");
362
363   Int_t det, t;
364   Float_t x, y, z, q, dy, dydx, dzdx, cov[3], covcl[3];
365   TH2I *h2 = 0x0;
366
367   // define limits around ExB for which x contribution is negligible
368   const Float_t kDtgPhi = 3.5e-2; //(+- 2 deg)
369
370   TObjArray *arr0 = (TObjArray*)fContainer->At(kCenter);
371   TObjArray *arr1 = (TObjArray*)fContainer->At(kSigm);
372   TObjArray *arr2 = (TObjArray*)fContainer->At(kMean);
373
374   const AliTRDclusterInfo *cli = 0x0;
375   TIterator *iter=fInfo->MakeIterator();
376   while((cli=dynamic_cast<AliTRDclusterInfo*>((*iter)()))){
377     cli->GetCluster(det, x, y, z, q, t, covcl);
378     if(fDet>=0 && fDet!=det) continue;
379     
380     dy = cli->GetResolution();
381     cli->GetGlobalPosition(y, z, dydx, dzdx, &cov[0]);
382
383     // resolution as a function of cluster charge
384     // only for phi equal exB 
385     if(TMath::Abs(dydx-fExB) < kDtgPhi){
386       h2 = (TH2I*)fContainer->At(kQRes);
387       h2->Fill(TMath::Log(q), dy);
388     }
389
390     // do not use problematic clusters in resolution analysis
391     // TODO define limits as calibration aware (gain) !!
392     if(q<20. || q>250.) continue;
393
394     // resolution as a function of y displacement from pad center
395     // only for phi equal exB
396     if(TMath::Abs(dydx-fExB) < kDtgPhi){
397       h2 = (TH2I*)arr0->At(AliTRDgeometry::GetLayer(det));
398       h2->Fill(cli->GetYDisplacement(), dy);
399     }
400
401     Int_t it = fAt->FindBin(cli->GetDriftLength());
402     if(it==0 || it == fAt->GetNbins()+1){
403       AliWarning(Form("Drift length %f outside allowed range", cli->GetDriftLength()));
404       continue;
405     }
406     Int_t id = fAd->FindBin(cli->GetAnisochronity());
407     if(id==0 || id == fAd->GetNbins()+1){
408       AliWarning(Form("Distance to anode %f outside allowed range", cli->GetAnisochronity()));
409       continue;
410     }
411     // calculate index of cluster position in array
412     Int_t hid = (id-1)*kNTB+it-1;
413
414     // fill histo for resolution (sigma)
415     h2 = (TH2I*)arr1->At(hid);
416     h2->Fill(dydx, dy);
417
418     // fill histo for systematic (mean)
419     h2 = (TH2I*)arr2->At(hid); 
420     h2->Fill(dydx-cli->GetTilt()*dzdx, dy);  
421   }
422   PostData(0, fContainer);
423 }
424
425
426 //_______________________________________________________
427 Bool_t AliTRDclusterResolution::PostProcess()
428 {
429   if(!fContainer) return kFALSE;
430   if(!HasExB()) AliWarning("ExB was not set. Call SetExB() before running the post processing.");
431   
432   TObjArray *arr = 0x0;
433   TH2 *h2 = 0x0; 
434   if(!fResults){
435     TGraphErrors *g = 0x0;
436     fResults = new TObjArray(sizeof(AliTRDclusterResolution::EResultContainers));
437     fResults->SetOwner();
438     fResults->AddAt(arr = new TObjArray(3), kQRes);
439     arr->SetOwner();
440     arr->AddAt(g = new TGraphErrors(), 0);
441     g->SetLineColor(kBlue); g->SetMarkerColor(kBlue);
442     g->SetMarkerStyle(7); 
443     arr->AddAt(g = new TGraphErrors(), 1);
444     g->SetLineColor(kRed); g->SetMarkerColor(kRed);
445     g->SetMarkerStyle(23); 
446     arr->AddAt(g = new TGraphErrors(), 2);
447     g->SetLineColor(kGreen); g->SetMarkerColor(kGreen);
448     g->SetMarkerStyle(7); 
449
450     fResults->AddAt(arr = new TObjArray(3), kCenter);
451     arr->SetOwner();
452
453     if(!(h2 = (TH2F*)gROOT->FindObject("hYM"))){
454       h2 = new TH2F("hYM", "", 
455       AliTRDgeometry::kNlayer, -.5, AliTRDgeometry::kNlayer-.5, 51, -.51, .51);
456     }
457     arr->AddAt(h2, 0);
458     h2->SetXTitle("ly");
459     h2->SetYTitle("y [w]");
460     h2->SetZTitle("#mu_{x} [cm]");
461     arr->AddAt(h2 = (TH2F*)h2->Clone("hYS"), 1);
462     h2->SetZTitle("#sigma_{x} [cm]");
463     arr->AddAt(h2 = (TH2F*)h2->Clone("hYP"), 2);
464     h2->SetZTitle("entries");
465
466     fResults->AddAt(arr = new TObjArray(2), kSigm);
467     arr->SetOwner();
468     if(!(h2 = (TH2F*)gROOT->FindObject("hSX"))){
469       h2 = new TH2F("hSX", "", 
470       fAd->GetNbins(), fAd->GetXmin(), fAd->GetXmax(), 
471       fAt->GetNbins(), fAt->GetXmin(), fAt->GetXmax());
472     }
473     arr->AddAt(h2, 0);
474     h2->SetXTitle("d [mm]");
475     h2->SetYTitle("x_{drift} [cm]");
476     h2->SetZTitle("#sigma_{x} [cm]");
477     arr->AddAt(h2 = (TH2F*)h2->Clone("hSY"), 1);
478     h2->SetZTitle("#sigma_{y} [cm]");
479
480     fResults->AddAt(arr = new TObjArray(2), kMean);
481     arr->SetOwner();
482     arr->AddAt(h2 = (TH2F*)h2->Clone("hDX"), 0);
483     h2->SetZTitle("dx [cm]");
484     arr->AddAt(h2 = (TH2F*)h2->Clone("hX0"), 1);
485     h2->SetZTitle("x0 [cm]");
486   } else {
487     TObject *o = 0x0;
488     TIterator *iter=fResults->MakeIterator();
489     while((o=((*iter)()))) o->Clear(); // maybe it is wrong but we should never reach this point
490   }
491   
492   // precalculated value of tg^2(alpha_L)
493   Double_t exb2 = fExB*fExB;
494   // square of the mean value of sigma drift length.
495   // has to come from previous calibration 
496   //Double_t sxd2 = 1.;// [mm^2]
497
498   printf("ExB[%e] ExB2[%e]\n", fExB, exb2);
499
500   // process resolution dependency on charge
501   if(HasProcessCharge()) ProcessCharge();
502   
503   // process resolution dependency on y displacement
504   if(HasProcessCenterPad()) ProcessCenterPad();
505
506   // process resolution dependency on drift legth and drift cell width
507   if(HasProcessSigma()) ProcessSigma();
508
509   // process systematic shift on drift legth and drift cell width
510   if(HasProcessMean()) ProcessMean();
511
512   return kTRUE;
513 }
514
515 //_______________________________________________________
516 Bool_t AliTRDclusterResolution::SetExB(Int_t det, Int_t col, Int_t row)
517 {
518   // check OCDB
519   AliCDBManager *cdb = AliCDBManager::Instance();
520   if(cdb->GetRun() < 0){
521     AliError("OCDB manager not properly initialized");
522     return kFALSE;
523   }
524
525   // check magnetic field
526   if(TMath::Abs(AliTracker::GetBz()) < 1.e-10){
527     AliWarning("B=0. Magnetic field may not be initialized. Continue if you know what you are doing !");
528   }
529
530   // set reference detector if any
531   if(det>=0 && det<AliTRDgeometry::kNdet) fDet = det;
532   else det = 0;
533
534   AliTRDcalibDB *fCalibration  = AliTRDcalibDB::Instance();
535   AliTRDCalROC  *fCalVdriftROC = fCalibration->GetVdriftROC(det);
536   const AliTRDCalDet  *fCalVdriftDet = fCalibration->GetVdriftDet();
537
538   fVdrift = fCalVdriftDet->GetValue(det) * fCalVdriftROC->GetValue(col, row);
539   fExB   = fCalibration->GetOmegaTau(fVdrift, -0.1*AliTracker::GetBz());
540   SetBit(kExB);
541   return kTRUE;
542 }
543
544 //_______________________________________________________
545 void AliTRDclusterResolution::SetVisual()
546 {
547   if(fCanvas) return;
548   fCanvas = new TCanvas("clResCanvas", "Cluster Resolution Visualization", 10, 10, 600, 600);
549 }
550
551 //_______________________________________________________
552 void AliTRDclusterResolution::ProcessCharge()
553 {
554   TH2I *h2 = 0x0;
555   if((h2 = (TH2I*)fContainer->At(kQRes))) {
556     AliWarning("Missing dy=f(Q) histo");
557     return;
558   }
559   TF1 f("f", "gaus", -.5, .5);
560   TAxis *ax = 0x0;
561   TH1D *h1 = 0x0;
562
563   TObjArray *arr = (TObjArray*)fResults->At(kQRes);
564   TGraphErrors *gqm = (TGraphErrors*)arr->At(0);
565   TGraphErrors *gqs = (TGraphErrors*)arr->At(1);
566   TGraphErrors *gqp = (TGraphErrors*)arr->At(2);
567   Double_t q, n = 0., entries;
568   ax = h2->GetXaxis();
569   for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
570     Double_t q = TMath::Exp(ax->GetBinCenter(ix));
571     if(q<20. || q>250.) continue; // ?!
572
573     h1 = h2->ProjectionY("py", ix, ix);
574     entries = h1->GetEntries();
575     if(entries < 50) continue;
576     Adjust(&f, h1);
577     h1->Fit(&f, "Q");
578
579     // Fill sy^2 = f(q)
580     Int_t ip = gqm->GetN();
581     gqm->SetPoint(ip, q, 10.*f.GetParameter(1));
582     gqm->SetPointError(ip, 0., 10.*f.GetParError(1));
583
584     // correct sigma for ExB effect
585     gqs->SetPoint(ip, q, 1.e1*f.GetParameter(2)/**f.GetParameter(2)-exb2*sxd2*/);
586     gqs->SetPointError(ip, 0., 1.e1*f.GetParError(2)/**f.GetParameter(2)*/);
587
588     // save probability
589     n += entries;
590     gqp->SetPoint(ip, q, entries);
591     gqp->SetPointError(ip, 0., 0./*TMath::Sqrt(entries)*/);
592   } 
593
594   // normalize probability and get mean sy
595   Double_t sm = 0., sy;
596   for(Int_t ip=gqp->GetN(); ip--;){
597     gqp->GetPoint(ip, q, entries);
598     entries/=n;
599     gqp->SetPoint(ip, q, entries);
600     gqs->GetPoint(ip, q, sy);
601     sm += entries*sy;
602   }
603
604   // error parametrization s(q) = <sy> + b(1/q-1/q0)
605   TF1 fq("fq", "[0] + [1]/x", 20., 250.);
606   gqs->Fit(&fq);
607   printf("sy(Q) :: sm[%f] b[%f] 1/q0[%f]\n", sm, fq.GetParameter(1), (sm-fq.GetParameter(0))/fq.GetParameter(1));
608 }
609
610 //_______________________________________________________
611 void AliTRDclusterResolution::ProcessCenterPad()
612 {
613   TObjArray *arr = (TObjArray*)fContainer->At(kCenter);
614   if(!arr) {
615     AliWarning("Missing dy=f(y_d) container");
616     return;
617   }
618   TF1 f("f", "gaus", -.5, .5);
619   TAxis *ax = 0x0;
620   TH1D *h1 = 0x0, *h11 = 0x0;
621   TH2I *h2 = 0x0;
622
623   TObjArray *arrg = (TObjArray*)fResults->At(kCenter);
624   TH2F *hym = (TH2F*)arrg->At(0);
625   TH2F *hys = (TH2F*)arrg->At(1);
626   TH2F *hyp = (TH2F*)arrg->At(2);
627   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
628     if(!(h2 = (TH2I*)arr->At(ily))) continue;
629     ax = h2->GetXaxis();
630     for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
631       Float_t yd = ax->GetBinCenter(ix);
632       h1 = h2->ProjectionY("py", ix, ix);
633       Int_t entries = (Int_t)h1->GetEntries();
634       if(entries < 50) continue;
635       //Adjust(&f, h1);
636       h1->Fit(&f, "QN");
637   
638       // Fill sy = f(y_w)
639       hyp->Fill(ily, yd, entries);
640       hym->Fill(ily, yd, f.GetParameter(1));
641       //hym->SetPointError(ip, 0., f.GetParError(1));
642       hys->Fill(ily, yd, f.GetParameter(2));
643       //hys->SetPointError(ip, 0., f.GetParError(2));
644     } 
645   }
646
647   // POSTPROCESS SPECTRA
648   // Found correction for systematic deviation  
649   TF1 fprf("fprf", "[0]+[1]*sin([2]*x)", -.5, .5);
650   fprf.SetParameter(0, 0.);
651   fprf.SetParameter(1, 1.1E-2);
652   fprf.SetParameter(2, -TMath::PiOver2()/0.25);
653   printf("  const Float_t cy[AliTRDgeometry::kNlayer][3] = {\n");
654   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
655     h1 = hym->ProjectionY("hym_pyy", ily+1, ily+1);
656     // adjust errors 
657     for(Int_t ib=h1->GetNbinsX(); ib--;) h1->SetBinError(ib, 0.002);
658     h1->Fit(&fprf, "Q");
659     printf("    {%5.3e, %5.3e, %5.3e},\n", fprf.GetParameter(0), fprf.GetParameter(1), fprf.GetParameter(2));
660
661     if(!fCanvas) continue;
662     fCanvas->cd(); 
663     h1->SetMinimum(-0.02);h1->SetMaximum(0.02);h1->Draw("e1"); 
664     h11 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
665     h11->Scale(.8/h11->Integral());
666     h11->SetLineColor(kBlue); h11->Draw("csame");
667     fCanvas->Modified(); fCanvas->Update(); 
668     if(IsSaveAs())
669     fCanvas->SaveAs(Form("Figures/ProcessCenterPad_M_Ly%d.gif", ily));
670     else gSystem->Sleep(500);
671   }
672   printf("  };\n");
673
674   // Parameterization for sigma PRF  
675   TF1 fgaus("fgaus", "gaus", -.5, .5);
676   printf("  const Float_t scy[AliTRDgeometry::kNlayer][4] = {\n");
677   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
678     h1 = hys->ProjectionY("hys_pyy", ily+1, ily+1);
679     // adjust errors 
680     for(Int_t ib=h1->GetNbinsX(); ib--;) h1->SetBinError(ib, 0.002);
681
682     h1->Fit(&fgaus, "Q");
683     printf("    {%5.3e, %5.3e, %5.3e, ", fgaus.GetParameter(0), fgaus.GetParameter(1), fgaus.GetParameter(2));
684
685     // calculate mean sigma on the pad center distribution
686     Float_t sy = 0.;
687     h1 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
688     for(Int_t ix=1; ix<=h1->GetNbinsX(); ix++){
689       sy += fgaus.Eval(h1->GetBinCenter(ix))*h1->GetBinContent(ix);
690     }
691     sy /= h1->GetEntries();
692     printf("%5.3e},\n", sy);
693
694     if(!fCanvas) continue;
695     fCanvas->cd(); 
696     h1->SetMinimum(0.01);h1->SetMaximum(0.04);h1->Draw("e1"); 
697     h11 = hyp->ProjectionY("hyp_pyy", ily+1, ily+1);
698     h11->Scale(1./h11->Integral());
699     h11->SetLineColor(kBlue); h11->Draw("csame");
700     fCanvas->Modified(); fCanvas->Update(); 
701     if(IsSaveAs())
702     fCanvas->SaveAs(Form("Figures/ProcessCenterPad_S_Ly%d.gif", ily));
703     else gSystem->Sleep(500);
704   }
705   printf("  };\n");
706 }
707
708 //_______________________________________________________
709 void AliTRDclusterResolution::ProcessSigma()
710 {
711   TObjArray *arr = (TObjArray*)fContainer->At(kSigm);
712   if(!arr){
713     AliWarning("Missing dy=f(x_d, d_wire) container");
714     return;
715   }
716   TLinearFitter gs(1,"pol1");
717   TF1 f("f", "gaus", -.5, .5);
718   TAxis *ax = 0x0;
719   TH1D *h1 = 0x0;
720   TH2I *h2 = 0x0;
721   // init visualization
722   TGraphErrors *ggs = 0x0;
723   TLine *line = 0x0;
724   if(fCanvas){
725     ggs = new TGraphErrors();
726     line = new TLine();
727   }
728
729   Double_t d(0.), x(0.), sx, sy;
730   TObjArray *arrr = (TObjArray*)fResults->At(kSigm);
731   TH2F *hsx = (TH2F*)arrr->At(0);
732   TH2F *hsy = (TH2F*)arrr->At(1);
733   for(Int_t id=1; id<=fAd->GetNbins(); id++){
734     d = fAd->GetBinCenter(id); //[mm]
735     printf(" Doing d = %5.3f [mm]\n", d);
736     for(Int_t it=1; it<=fAt->GetNbins(); it++){
737       x = fAt->GetBinCenter(it); //[mm]
738       Int_t idx = (id-1)*kNTB+it-1;
739
740       // retrieve data histogram
741       if(!(h2 = (TH2I*)arr->At(idx))) {
742         AliWarning(Form("Missing histo at index idx[%3d] [id[%2d] it[%2d]] xd[%f] d[%f]\n", idx, id, it, x, d));
743         continue;
744       }
745
746       if(fCanvas){ 
747         new(ggs) TGraphErrors();
748         ggs->SetMarkerStyle(7);
749       }
750       gs.ClearPoints();
751       ax = h2->GetXaxis();
752       for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
753         Float_t dydx = ax->GetBinCenter(ix);
754         Double_t tgg = (dydx-fExB)/(1.+dydx*fExB);
755         Double_t tgg2 = tgg*tgg;
756         h1 = h2->ProjectionY("py", ix, ix);
757         if(h1->GetEntries() < 100) continue;
758         //Adjust(&f, h1);
759         //printf("\tFit ix[%d] on %s entries [%d]\n", ix, h2->GetName(), (Int_t)h1->GetEntries());
760         h1->Fit(&f, "QN");
761         Double_t s2  = f.GetParameter(2)*f.GetParameter(2);
762         Double_t s2e = 2.*f.GetParameter(2)*f.GetParError(2);
763         // Fill sy^2 = f(tg^2(phi-a_L))
764         gs.AddPoint(&tgg2, s2, s2e);
765
766         if(!ggs) continue;
767         Int_t ip = ggs->GetN();
768         ggs->SetPoint(ip, tgg2, s2);
769         ggs->SetPointError(ip, 0., s2e);
770         //printf("%f, %f, %f, \n", tgg2, s2, s2e);
771       }
772       if(gs.Eval()) continue;
773       sx = gs.GetParameter(1); 
774       sy = gs.GetParameter(0);
775       hsx->SetBinContent(id, it, TMath::Sqrt(sx));
776       //hsx->SetBinError(id, it, sig.GetParError(1));
777
778       // s^2_y = s0^2_y + tg^2(a_L) * s^2_x 
779       hsy->SetBinContent(id, it, TMath::Sqrt(sy)/* - exb2*sig.GetParameter(1)*/);
780       //hsy->SetBinError(id, it, sig.GetParError(0)+exb2*exb2*sig.GetParError(1));
781
782       /*if(!fCanvas) */continue;
783       fCanvas->cd();
784       new(line) TLine(0, sy, .025, sy+.025*sx); 
785       line->SetLineColor(kRed);line->SetLineWidth(2);
786       ggs->Draw("apl"); line->Draw("");
787       fCanvas->Modified(); fCanvas->Update();
788       if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessSigma_D%d_T%02d.gif", id, it));
789       else gSystem->Sleep(100);
790
791       printf("    xd=%4.1f[cm] sx=%5.3e[cm] sy=%5.3e[cm]\n", x, sx, sy);
792     }
793   }
794 //   if(ggs) delete ggs; // TODO fix this memory leak !
795 //   if(line) delete line;
796
797   // fit sy parallel to the drift
798   TF1 fsyD("fsy", "[0]+[1]*exp(1./(x+[2]))", 0.5, 3.5);
799   fsyD.SetParameter(0, .025);
800   fsyD.SetParameter(1, -8.e-4);
801   fsyD.SetParameter(2, -.25);
802   h1 = hsy->ProjectionY("hsy_pyy"); h1->Scale(1./hsy->GetNbinsX());
803   for(Int_t ib=1; ib<=h1->GetNbinsX(); ib++) h1->SetBinError(ib, 0.002);
804   h1->Fit(&fsyD, "Q", "", .5, 3.5);
805   printf("  const Float_t sy0 = %5.3e;\n", fsyD.GetParameter(0));
806   printf("  const Float_t sya = %5.3e;\n", fsyD.GetParameter(1));
807   printf("  const Float_t syb = %5.3e;\n", fsyD.GetParameter(2));
808   if(fCanvas){
809     fCanvas->cd();
810     h1->Draw("e1"); fsyD.Draw("same");
811     fCanvas->Modified(); fCanvas->Update();
812     if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SY||.gif");
813     else gSystem->Sleep(1000);
814   }
815   
816
817   // fit sx parallel to the drift
818   TF1 fsxD("fsx", "gaus(0)+expo(3)", 0., 3.5);
819   h1 = hsx->ProjectionY("hsx_pyy"); 
820   h1->Scale(1./hsx->GetNbinsX());
821   for(Int_t ib=1; ib<=h1->GetNbinsX(); ib++) if(h1->GetBinContent(ib)>0.) h1->SetBinError(ib, 0.02);
822   h1->Fit(&f, "QN", "", 0., 1.3);
823   fsxD.SetParameter(0, 0.);
824   fsxD.SetParameter(1, f.GetParameter(1));
825   fsxD.SetParameter(2, f.GetParameter(2));
826   TF1 expo("fexpo", "expo", 0., 3.5);
827   h1->Fit(&expo, "QN");
828   fsxD.SetParameter(3, expo.GetParameter(0));
829   fsxD.SetParameter(4, expo.GetParameter(1));
830   h1->Fit(&fsxD, "Q");
831   printf("  const Float_t sxgc = %5.3e;\n", fsxD.GetParameter(0));
832   printf("  const Float_t sxgm = %5.3e;\n", fsxD.GetParameter(1));
833   printf("  const Float_t sxgs = %5.3e;\n", fsxD.GetParameter(2));
834   printf("  const Float_t sxe0 = %5.3e;\n", fsxD.GetParameter(3));
835   printf("  const Float_t sxe1 = %5.3e;\n", fsxD.GetParameter(4));
836   if(fCanvas){
837     fCanvas->cd();
838     h1->Draw("e1"); fsxD.Draw("same");
839     fCanvas->Modified(); fCanvas->Update();
840     if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SX||.gif");
841     else gSystem->Sleep(1000);
842   }
843   
844   // fit sx perpendicular to the drift
845   Int_t nb = 0;
846   TAxis *ay = hsx->GetYaxis();
847   TH1D *hx = hsx->ProjectionX("hsx_px", 1,1);
848   TH1D *hAn = (TH1D*)hx->Clone("hAn"); hAn->Reset();
849   for(Int_t ib = 11; ib<24; ib++, nb++){
850     hx = hsx->ProjectionX("hsx_px", ib, ib);
851     for(Int_t id=1; id<=hx->GetNbinsX(); id++)
852       hx->SetBinContent(id, hx->GetBinContent(id)-fsxD.Eval(ay->GetBinCenter(ib))); 
853     hAn->Add(hx);
854   }
855   TF1 fsxA("fsxA", "pol2", 0., 0.25);
856   hAn->Scale(1./nb);
857   for(Int_t ib=1; ib<=hAn->GetNbinsX(); ib++) hAn->SetBinError(ib, 0.002);
858   hAn->Fit(&fsxA, "Q");
859   printf("  const Float_t sxd0 = %5.3e;\n", fsxA.GetParameter(0));
860   printf("  const Float_t sxd1 = %5.3e;\n", fsxA.GetParameter(1));
861   printf("  const Float_t sxd2 = %5.3e;\n", fsxA.GetParameter(2));
862   if(fCanvas){
863     fCanvas->cd();
864     hAn->Draw("e1"); fsxA.Draw("same");
865     fCanvas->Modified(); fCanvas->Update();
866     if(IsSaveAs()) fCanvas->SaveAs("Figures/ProcessSigma_SXL.gif");
867     else gSystem->Sleep(100);
868   }
869 }
870
871 //_______________________________________________________
872 void AliTRDclusterResolution::ProcessMean()
873 {
874   TObjArray *arr = (TObjArray*)fContainer->At(kMean);
875   if(!arr){
876     AliWarning("Missing dy=f(x_d, d_wire) container");
877     return;
878   }
879   TGraphErrors *gm = new TGraphErrors();
880   TF1 f("f", "gaus", -.5, .5);
881   TF1 line("l", Form("%6.4f*[0]+[1]*x", fExB), -.15, .15);
882   Double_t d(0.), x(0.), dx, x0;
883   TAxis *ax = 0x0;
884   TH1D *h1 = 0x0;
885   TH2I *h2 = 0x0;
886
887   TObjArray *arrr = (TObjArray*)fResults->At(kMean);
888   TH2F *hdx = (TH2F*)arrr->At(0);
889   TH2F *hx0 = (TH2F*)arrr->At(1);
890   for(Int_t id=1; id<=fAd->GetNbins(); id++){
891     d = fAd->GetBinCenter(id); //[mm]
892     printf(" Doing d = %5.3f [mm]\n", d);
893     for(Int_t it=1; it<=fAt->GetNbins(); it++){
894       x = fAt->GetBinCenter(it); //[mm]
895       Int_t idx = (id-1)*kNTB+it-1;
896       if(!(h2 = (TH2I*)arr->At(idx))) {
897         AliWarning(Form("Missing histo at index idx[%3d] [id[%2d] it[%2d]] xd[%f] d[%f]\n", idx, id, it, x, d));
898         continue;
899       }
900
901       new(gm) TGraphErrors();
902       gm->SetMarkerStyle(7);
903       ax = h2->GetXaxis();
904       for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
905         Double_t dydx = ax->GetBinCenter(ix);
906         h1 = h2->ProjectionY("py", ix, ix);
907         if(h1->GetEntries() < 50) continue;
908         //Adjust(&f, h1);
909         h1->Fit(&f, "QN");
910
911         // Fill <Dy> = f(dydx - h*dzdx)
912         Int_t ip = gm->GetN();
913         gm->SetPoint(ip, dydx, f.GetParameter(1));
914         gm->SetPointError(ip, 0., f.GetParError(1));
915         ip++;
916       }
917       if(gm->GetN()<4) continue;
918
919       gm->Fit(&line, "QN");
920       dx = line.GetParameter(1); // = (fVdrift + dVd)*(fT + tCorr);
921       x0 = line.GetParameter(0); // = tg(a_L)*dx
922       hdx->SetBinContent(id, it, dx);
923       hx0->SetBinContent(id, it, x0);
924
925       if(!fCanvas) continue;
926       fCanvas->cd();
927       gm->Draw("apl"); line.Draw("same");
928       fCanvas->Modified(); fCanvas->Update();
929       if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessMean_D%d_T%02d.gif", id, it));
930       else gSystem->Sleep(100);
931       printf("    xd=%4.1f[cm] dx=%5.3e[cm] x0=%5.3e[cm]\n", x, dx, x0);
932     }
933   }
934
935   TH1D *h=hdx->ProjectionY("hdx_py"); h->Scale(1./hdx->GetNbinsX());
936   printf("  const Float_t cx[] = {\n    ");
937   for(Int_t ib=1; ib<=h->GetNbinsX(); ib++){
938     printf("%6.4e, ", h->GetBinContent(ib));
939     if(ib%6==0) printf("\n    ");
940   }
941   printf("  };\n");
942
943   // delete gm; TODO memory leak ?
944 }