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