]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDtrackingResolution.cxx
fix PID reference figures style (AlexW)
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDtrackingResolution.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: AliTRDtrackingResolution.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  TRD tracking resolution                                               //
21 //
22 // The class performs resolution and residual studies 
23 // of the TRD tracks for the following quantities :
24 //   - spatial position (y, [z])
25 //   - angular (phi) tracklet
26 //   - momentum at the track level
27 // 
28 // The class has to be used for regular detector performance checks using the official macros:
29 //   - $ALICE_ROOT/TRD/qaRec/run.C
30 //   - $ALICE_ROOT/TRD/qaRec/makeResults.C
31 // 
32 // For stand alone usage please refer to the following example: 
33 // {  
34 //   gSystem->Load("libANALYSIS.so");
35 //   gSystem->Load("libTRDqaRec.so");
36 //   AliTRDtrackingResolution *res = new AliTRDtrackingResolution();
37 //   //res->SetMCdata();
38 //   //res->SetVerbose();
39 //   //res->SetVisual();
40 //   res->Load("TRD.TaskResolution.root");
41 //   if(!res->PostProcess()) return;
42 //   res->GetRefFigure(0);
43 // }  
44 //
45 //  Authors:                                                              //
46 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
47 //    Markus Fasel <M.Fasel@gsi.de>                                       //
48 //                                                                        //
49 ////////////////////////////////////////////////////////////////////////////
50
51 #include <cstring>
52
53 #include <TROOT.h>
54 #include <TSystem.h>
55 #include <TPDGCode.h>
56 #include <TObjArray.h>
57 #include <TH2.h>
58 #include <TH1.h>
59 #include <TF1.h>
60 #include <TCanvas.h>
61 #include <TBox.h>
62 #include <TProfile.h>
63 #include <TGraphErrors.h>
64 #include <TMath.h>
65 #include "TTreeStream.h"
66 #include "TGeoManager.h"
67
68 #include "AliAnalysisManager.h"
69 #include "AliTrackReference.h"
70 #include "AliTrackPointArray.h"
71 #include "AliCDBManager.h"
72
73 #include "AliTRDSimParam.h"
74 #include "AliTRDgeometry.h"
75 #include "AliTRDpadPlane.h"
76 #include "AliTRDcluster.h"
77 #include "AliTRDseedV1.h"
78 #include "AliTRDtrackV1.h"
79 #include "AliTRDtrackerV1.h"
80 #include "AliTRDReconstructor.h"
81 #include "AliTRDrecoParam.h"
82
83 #include "AliTRDtrackInfo/AliTRDclusterInfo.h"
84 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
85 #include "AliTRDtrackingResolution.h"
86
87 ClassImp(AliTRDtrackingResolution)
88
89 //________________________________________________________
90 AliTRDtrackingResolution::AliTRDtrackingResolution()
91   :AliTRDrecoTask("Resolution", "Tracking Resolution")
92   ,fStatus(0)
93   ,fReconstructor(0x0)
94   ,fGeo(0x0)
95   ,fGraphS(0x0)
96   ,fGraphM(0x0)
97   ,fCl(0x0)
98   ,fTrklt(0x0)
99   ,fMCcl(0x0)
100   ,fMCtrklt(0x0)
101 {
102   fReconstructor = new AliTRDReconstructor();
103   fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
104   fGeo = new AliTRDgeometry();
105
106   InitFunctorList();
107
108   DefineOutput(1, TObjArray::Class()); // cluster2track
109   DefineOutput(2, TObjArray::Class()); // tracklet2track
110   DefineOutput(3, TObjArray::Class()); // cluster2mc
111   DefineOutput(4, TObjArray::Class()); // tracklet2mc
112 }
113
114 //________________________________________________________
115 AliTRDtrackingResolution::~AliTRDtrackingResolution()
116 {
117   if(fGraphS){fGraphS->Delete(); delete fGraphS;}
118   if(fGraphM){fGraphM->Delete(); delete fGraphM;}
119   delete fGeo;
120   delete fReconstructor;
121   if(gGeoManager) delete gGeoManager;
122   if(fCl){fCl->Delete(); delete fCl;}
123   if(fTrklt){fTrklt->Delete(); delete fTrklt;}
124   if(fMCcl){fMCcl->Delete(); delete fMCcl;}
125   if(fMCtrklt){fMCtrklt->Delete(); delete fMCtrklt;}
126 }
127
128
129 //________________________________________________________
130 void AliTRDtrackingResolution::CreateOutputObjects()
131 {
132   // spatial resolution
133   OpenFile(0, "RECREATE");
134
135   fContainer = Histos();
136
137   fCl = new TObjArray();
138   fCl->SetOwner(kTRUE);
139   fTrklt = new TObjArray();
140   fTrklt->SetOwner(kTRUE);
141   fMCcl = new TObjArray();
142   fMCcl->SetOwner(kTRUE);
143   fMCtrklt = new TObjArray();
144   fMCtrklt->SetOwner(kTRUE);
145 }
146
147 //________________________________________________________
148 void AliTRDtrackingResolution::Exec(Option_t *opt)
149 {
150   fCl->Delete();
151   fTrklt->Delete();
152   fMCcl->Delete();
153   fMCtrklt->Delete();
154
155   AliTRDrecoTask::Exec(opt);
156
157   PostData(1, fCl);
158   PostData(2, fTrklt);
159   PostData(3, fMCcl);
160   PostData(4, fMCtrklt);
161 }
162
163 //________________________________________________________
164 TH1* AliTRDtrackingResolution::PlotCluster(const AliTRDtrackV1 *track)
165 {
166   if(track) fTrack = track;
167   if(!fTrack){
168     AliWarning("No Track defined.");
169     return 0x0;
170   }
171   TH1 *h = 0x0;
172   if(!(h = ((TH2I*)fContainer->At(kCluster)))){
173     AliWarning("No output histogram defined.");
174     return 0x0;
175   }
176
177   Double_t cov[3];
178   Float_t x0, y0, z0, dy, dydx, dzdx;
179   AliTRDseedV1 *fTracklet = 0x0;  
180   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
181     if(!(fTracklet = fTrack->GetTracklet(ily))) continue;
182     if(!fTracklet->IsOK()) continue;
183     x0 = fTracklet->GetX0();
184
185     // retrive the track angle with the chamber
186     y0   = fTracklet->GetYref(0);
187     z0   = fTracklet->GetZref(0);
188     dydx = fTracklet->GetYref(1);
189     dzdx = fTracklet->GetZref(1);
190     fTracklet->GetCovRef(cov);
191     Float_t tilt = fTracklet->GetTilt();
192     AliTRDcluster *c = 0x0;
193     fTracklet->ResetClusterIter(kFALSE);
194     while((c = fTracklet->PrevCluster())){
195       Float_t xc = c->GetX();
196       Float_t yc = c->GetY();
197       Float_t zc = c->GetZ();
198       Float_t dx = x0 - xc; 
199       Float_t yt = y0 - dx*dydx;
200       Float_t zt = z0 - dx*dzdx; 
201       yc -= tilt*(zc-zt); // tilt correction
202       dy = yt - yc;
203
204       h->Fill(dydx, dy/TMath::Sqrt(cov[0] + c->GetSigmaY2()));
205   
206       if(fDebugLevel>=1){
207         // Get z-position with respect to anode wire
208         AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
209         Int_t istk = fGeo->GetStack(c->GetDetector());
210         AliTRDpadPlane *pp = fGeo->GetPadPlane(ily, istk);
211         Float_t row0 = pp->GetRow0();
212         Float_t d  =  row0 - zt + pp->GetAnodeWireOffset();
213         d -= ((Int_t)(2 * d)) / 2.0;
214         if (d > 0.25) d  = 0.5 - d;
215
216 /*        AliTRDclusterInfo *clInfo = new AliTRDclusterInfo;
217         fCl->Add(clInfo);
218         clInfo->SetCluster(c);
219         clInfo->SetGlobalPosition(yt, zt, dydx, dzdx);
220         clInfo->SetResolution(dy);
221         clInfo->SetAnisochronity(d);
222         clInfo->SetDriftLength(dx);
223         (*fDebugStream) << "ClusterResiduals"
224           <<"clInfo.=" << clInfo
225           << "\n";*/
226       }
227     }
228   }
229   return h;
230 }
231
232
233 //________________________________________________________
234 TH1* AliTRDtrackingResolution::PlotTracklet(const AliTRDtrackV1 *track)
235 {
236   if(track) fTrack = track;
237   if(!fTrack){
238     AliWarning("No Track defined.");
239     return 0x0;
240   }
241   TH1 *h = 0x0;
242   if(!(h = ((TH2I*)fContainer->At(kTrackletY)))){
243     AliWarning("No output histogram defined.");
244     return 0x0;
245   }
246
247   Double_t cov[3], covR[3];
248   Float_t xref, y0, dx, dy, dydx;
249   AliTRDseedV1 *fTracklet = 0x0;  
250   for(Int_t il=AliTRDgeometry::kNlayer; il--;){
251     if(!(fTracklet = fTrack->GetTracklet(il))) continue;
252     if(!fTracklet->IsOK()) continue;
253     y0   = fTracklet->GetYref(0);
254     dydx = fTracklet->GetYref(1);
255     xref = fTracklet->GetXref();
256     dx   = fTracklet->GetX0() - xref;
257     dy   = y0-dx*dydx - fTracklet->GetYat(xref);
258     fTracklet->GetCovAt(xref, cov);
259     fTracklet->GetCovRef(covR);
260     h->Fill(dydx, dy/TMath::Sqrt(cov[0] + covR[0]));
261   }
262   return h;
263 }
264
265 //________________________________________________________
266 TH1* AliTRDtrackingResolution::PlotTrackletPhi(const AliTRDtrackV1 *track)
267 {
268   if(track) fTrack = track;
269   if(!fTrack){
270     AliWarning("No Track defined.");
271     return 0x0;
272   }
273   TH1 *h = 0x0;
274   if(!(h = ((TH2I*)fContainer->At(kTrackletPhi)))){
275     AliWarning("No output histogram defined.");
276     return 0x0;
277   }
278
279   AliTRDseedV1 *tracklet = 0x0;
280   for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++){
281     if(!(tracklet = fTrack->GetTracklet(il))) continue;
282     if(!tracklet->IsOK()) continue;
283     h->Fill(tracklet->GetYref(1), TMath::ATan(tracklet->GetYref(1))-TMath::ATan(tracklet->GetYfit(1)));
284   }
285   return h;
286 }
287
288
289 //________________________________________________________
290 TH1* AliTRDtrackingResolution::PlotResolution(const AliTRDtrackV1 *track)
291 {
292   if(!HasMCdata()){ 
293     AliWarning("No MC defined. Results will not be available.");
294     return 0x0;
295   }
296   if(track) fTrack = track;
297   if(!fTrack){
298     AliWarning("No Track defined.");
299     return 0x0;
300   }
301   TH1 *h = 0x0;
302   UChar_t s;
303   Int_t pdg = fMC->GetPDG(), det=-1;
304   Int_t label = fMC->GetLabel();
305   Float_t p, pt, xref, x0, y0, z0, dx, dy, dz, dydx, dzdx;
306
307   if(fDebugLevel>=1){
308     Double_t DX[12], DY[12], DZ[12], DPt[12], COV[12][15];
309     fMC->PropagateKalman(DX, DY, DZ, DPt, COV);
310     (*fDebugStream) << "MCkalman"
311       << "pdg="  << pdg
312       << "dx0="  << DX[0]
313       << "dx1="  << DX[1]
314       << "dx2="  << DX[2]
315       << "dy0="  << DY[0]
316       << "dy1="  << DY[1]
317       << "dy2="  << DY[2]
318       << "dz0="  << DZ[0]
319       << "dz1="  << DZ[1]
320       << "dz2="  << DZ[2]
321       << "dpt0=" << DPt[0]
322       << "dpt1=" << DPt[1]
323       << "dpt2=" << DPt[2]
324       << "\n";
325   }
326
327   AliTRDseedV1 *fTracklet = 0x0;  
328   for(Int_t ily=0; ily<AliTRDgeometry::kNlayer; ily++){
329     if(!(fTracklet = fTrack->GetTracklet(ily)) ||
330        !fTracklet->IsOK()) continue;
331
332     det = fTracklet->GetDetector();
333     x0  = fTracklet->GetX0();
334     //radial shift with respect to the MC reference (radial position of the pad plane)
335     xref= fTracklet->GetXref();
336     dx  = x0 - xref;
337     if(!fMC->GetDirections(x0, y0, z0, dydx, dzdx, pt, s)) continue;
338     // MC track position at reference radial position
339     Float_t yt = y0 - (x0-xref)*dydx;
340     Float_t zt = z0 - (x0-xref)*dzdx;
341     p = pt*(1.+dzdx*dzdx); // pt -> p
342
343     // add Kalman residuals for y, z and pt
344     Float_t yr = fTracklet->GetYref(0) - dx*fTracklet->GetYref(1);
345     dy = yt - yr;
346     Float_t zr = fTracklet->GetZref(0) - dx*fTracklet->GetZref(1);
347     dz = zt - zr;
348     Float_t tgl = fTracklet->GetTgl();
349     Float_t ptr = fTracklet->GetMomentum()/(1.+tgl*tgl);
350
351     ((TH2I*)fContainer->At(kMCtrackY))->Fill(dydx, dy);
352     ((TH2I*)fContainer->At(kMCtrackZ))->Fill(dzdx, dz);
353     if(pdg!=kElectron && pdg!=kPositron) ((TH2I*)fContainer->At(kMCtrackPt))->Fill(1./pt, ptr-pt);
354     // Fill Debug stream for Kalman track
355     if(fDebugLevel>=1){
356       Float_t dydxr = fTracklet->GetYref(1);
357       (*fDebugStream) << "MCtrack"
358         << "det="     << det
359         << "pdg="     << pdg
360         << "pt="      << pt
361         << "yt="      << yt
362         << "zt="      << zt
363         << "dydx="    << dydx
364         << "dzdx="    << dzdx
365         << "ptr="     << ptr
366         << "dy="      << dy
367         << "dz="      << dz
368         << "dydxr="   << dydxr
369         << "dzdxr="   << tgl
370         << "\n";
371     }
372
373     // recalculate tracklet based on the MC info
374     AliTRDseedV1 tt(*fTracklet);
375     tt.SetZref(0, z0);
376     tt.SetZref(1, dzdx); 
377     if(!tt.Fit(kTRUE)) continue;
378     xref= fTracklet->GetXref(); // the true one 
379     yt = y0 - (x0-xref)*dydx;
380
381     // add tracklet residuals for y and dydx
382     Float_t yf = tt.GetYat(xref);
383     dy = yt - yf;
384     Float_t dphi   = (tt.GetYfit(1) - dydx);
385     dphi /= 1.- tt.GetYfit(1)*dydx;
386     ((TH2I*)fContainer->At(kMCtrackletY))->Fill(dydx, dy);
387     ((TH2I*)fContainer->At(kMCtrackletPhi))->Fill(dydx, dphi*TMath::RadToDeg());
388
389     Float_t dz = 100.;
390     Bool_t rc = fTracklet->IsRowCross(); 
391     if(rc){
392       // add tracklet residuals for z
393       Double_t *xyz = tt.GetCrossXYZ();
394       dz = xyz[2] - (z0 - (x0 - xyz[0])*dzdx) ;
395       ((TH2I*)fContainer->At(kMCtrackletZ))->Fill(dzdx, dz);
396     }
397   
398     // Fill Debug stream for tracklet
399     if(fDebugLevel>=1){
400       (*fDebugStream) << "MCtracklet"
401         << "det="     << det
402         << "pdg="     << pdg
403         << "p="       << p
404         << "yt="      << yt
405         << "zt="      << zt
406         << "dydx="    << dydx
407         << "dzdx="    << dzdx
408         << "rowc="    << rc
409         << "dy="      << dy
410         << "dz="      << dz
411         << "dphi="    << dphi
412         << "\n";
413     }
414
415     Int_t istk = AliTRDgeometry::GetStack(det); 
416     AliTRDpadPlane *pp = fGeo->GetPadPlane(ily, istk);
417     Float_t zr0 = pp->GetRow0() + pp->GetAnodeWireOffset();
418     Float_t tilt = fTracklet->GetTilt();
419
420     Double_t x,y;
421     AliTRDcluster *c = 0x0;
422     fTracklet->ResetClusterIter(kFALSE);
423     while((c = fTracklet->PrevCluster())){
424       Float_t  q = TMath::Abs(c->GetQ());
425       //AliTRDseedV1::GetClusterXY(c,x,y);
426       x = c->GetX(); y = c->GetY();
427       Float_t xc = x;
428       Float_t yc = y;
429       Float_t zc = c->GetZ();
430       dx = x0 - xc; 
431       yt = y0 - dx*dydx;
432       zt = z0 - dx*dzdx;
433       dy = yt - (yc - tilt*(zc-zt));
434
435       // Fill Histograms
436       if(q>20. && q<250.) ((TH2I*)fContainer->At(kMCcluster))->Fill(dydx, dy);
437
438       // Fill calibration container
439       Float_t d = zr0 - zt;
440       d -= ((Int_t)(2 * d)) / 2.0;
441       if (d > 0.25) d  = 0.5 - d;
442       AliTRDclusterInfo *clInfo = new AliTRDclusterInfo;
443       fMCcl->Add(clInfo);
444       clInfo->SetCluster(c);
445       clInfo->SetMC(pdg, label);
446       clInfo->SetGlobalPosition(yt, zt, dydx, dzdx);
447       clInfo->SetResolution(dy);
448       clInfo->SetAnisochronity(d);
449       clInfo->SetDriftLength(dx-.5*AliTRDgeometry::CamHght());
450       clInfo->SetTilt(tilt);
451
452       // Fill Debug Tree
453       if(fDebugLevel>=2){
454         //clInfo->Print();
455         (*fDebugStream) << "MCcluster"
456           <<"clInfo.=" << clInfo
457           << "\n";
458       }
459     }
460   }
461   return h;
462 }
463
464
465 //________________________________________________________
466 Bool_t AliTRDtrackingResolution::GetRefFigure(Int_t ifig)
467 {
468   Float_t y[2] = {0., 0.};
469   TBox *b = 0x0;
470   TAxis *ax = 0x0;
471   TGraphErrors *g = 0x0;
472   switch(ifig){
473   case kCluster:
474     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
475     g->Draw("apl");
476     ax = g->GetHistogram()->GetYaxis();
477     y[0] = -0.5; y[1] = 2.5;
478     ax->SetRangeUser(y[0], y[1]);
479     ax->SetTitle("Cluster-Track Pools #sigma/#mu [mm]");
480     ax = g->GetHistogram()->GetXaxis();
481     ax->SetTitle("tg(#phi)");
482     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
483     g->Draw("pl");
484     b = new TBox(-.15, y[0], .15, y[1]);
485     b->SetFillStyle(3002);b->SetFillColor(kGreen);
486     b->SetLineColor(0); b->Draw();
487     return kTRUE;
488   case kTrackletY:
489     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
490     g->Draw("apl");
491     ax = g->GetHistogram()->GetYaxis();
492     ax->SetRangeUser(-.5, 3.);
493     ax->SetTitle("Tracklet-Track Y-Pools #sigma/#mu [mm]");
494     ax = g->GetHistogram()->GetXaxis();
495     ax->SetTitle("tg(#phi)");
496     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
497     g->Draw("pl");
498     b = new TBox(-.15, y[0], .15, y[1]);
499     b->SetFillStyle(3002);b->SetFillColor(kGreen);
500     b->SetLineColor(0); b->Draw();
501     return kTRUE;
502   case kTrackletPhi:
503     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
504     g->Draw("apl");
505     ax = g->GetHistogram()->GetYaxis();
506     ax->SetRangeUser(-.5, 2.);
507     ax->SetTitle("Tracklet-Track Phi-Residuals #sigma/#mu [rad]");
508     ax = g->GetHistogram()->GetXaxis();
509     ax->SetTitle("tg(#phi)");
510     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
511     g->Draw("pl");
512     b = new TBox(-.15, y[0], .15, y[1]);
513     b->SetFillStyle(3002);b->SetFillColor(kGreen);
514     b->SetLineColor(0); b->Draw();
515     return kTRUE;
516   case kMCcluster:
517     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
518     ax = g->GetHistogram()->GetYaxis();
519     y[0] = -.1; y[1] = 1.5;
520     ax->SetRangeUser(y[0], y[1]);
521     ax->SetTitle("Y_{cluster} #sigma/#mu [mm]");
522     ax = g->GetHistogram()->GetXaxis();
523     ax->SetTitle("tg(#phi)");
524     g->Draw("apl");
525     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
526     g->Draw("pl");
527     b = new TBox(-.15, y[0], .15, y[1]);
528     b->SetFillStyle(3002);b->SetFillColor(kBlue);
529     b->SetLineColor(0); b->Draw();
530     return kTRUE;
531   case kMCtrackletY:
532     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
533     ax = g->GetHistogram()->GetYaxis();
534     y[0] = -.05; y[1] = 0.3;
535     ax->SetRangeUser(y[0], y[1]);
536     ax->SetTitle("Y_{tracklet} #sigma/#mu [mm]");
537     ax = g->GetHistogram()->GetXaxis();
538     ax->SetTitle("tg(#phi)");
539     g->Draw("apl");
540     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
541     g->Draw("pl");
542     b = new TBox(-.15, y[0], .15, y[1]);
543     b->SetFillStyle(3002);b->SetFillColor(kBlue);
544     b->SetLineColor(0); b->Draw();
545     return kTRUE;
546   case kMCtrackletZ:
547     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
548     ax = g->GetHistogram()->GetYaxis();
549     ax->SetRangeUser(-.5, 1.);
550     ax->SetTitle("Z_{tracklet} #sigma/#mu [mm]");
551     ax = g->GetHistogram()->GetXaxis();
552     ax->SetTitle("tg(#theta)");
553     g->Draw("apl");
554     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
555     g->Draw("pl");
556     return kTRUE;
557   case kMCtrackletPhi:
558     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
559     ax = g->GetHistogram()->GetYaxis();
560     y[0] = -.05; y[1] = .2;
561     ax->SetRangeUser(y[0], y[1]);
562     ax->SetTitle("#Phi_{tracklet} #sigma/#mu [deg]");
563     ax = g->GetHistogram()->GetXaxis();
564     ax->SetTitle("tg(#phi)");
565     g->Draw("apl");
566     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
567     g->Draw("pl");
568     return kTRUE;
569   case kMCtrackY:
570     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
571     ax = g->GetHistogram()->GetYaxis();
572     y[0] = -.05; y[1] = 0.2;
573     ax->SetRangeUser(y[0], y[1]);
574     ax->SetTitle("Y_{track} #sigma/#mu [mm]");
575     ax = g->GetHistogram()->GetXaxis();
576     ax->SetTitle("tg(#phi)");
577     g->Draw("apl");
578     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
579     g->Draw("pl");
580     b = new TBox(-.15, y[0], .15, y[1]);
581     b->SetFillStyle(3002);b->SetFillColor(kBlue);
582     b->SetLineColor(0); b->Draw();
583     return kTRUE;
584   case kMCtrackZ:
585     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
586     ax = g->GetHistogram()->GetYaxis();
587     ax->SetRangeUser(-.5, 2.);
588     ax->SetTitle("Z_{track} #sigma/#mu [mm]");
589     ax = g->GetHistogram()->GetXaxis();
590     ax->SetTitle("tg(#theta)");
591     g->Draw("apl");
592     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
593     g->Draw("pl");
594     return kTRUE;
595   case kMCtrackPt:
596     if(!(g = (TGraphErrors*)fGraphS->At(ifig))) break;
597     ax = g->GetHistogram()->GetYaxis();
598     ax->SetRangeUser(-.5, 2.);
599     ax->SetTitle("#epsilon_{P_{t}}^{track} / #mu [%]");
600     ax = g->GetHistogram()->GetXaxis();
601     ax->SetTitle("1/p_{t}");
602     g->Draw("apl");
603     if(!(g = (TGraphErrors*)fGraphM->At(ifig))) break;
604     g->Draw("pl");
605     return kTRUE;
606   }
607   AliInfo(Form("Reference plot [%d] missing result", ifig));
608   return kFALSE;
609 }
610
611
612 //________________________________________________________
613 Bool_t AliTRDtrackingResolution::PostProcess()
614 {
615   //fContainer = dynamic_cast<TObjArray*>(GetOutputData(0));
616   if (!fContainer) {
617     Printf("ERROR: list not available");
618     return kFALSE;
619   }
620   fNRefFigures = fContainer->GetEntriesFast();
621   TGraphErrors *gm = 0x0, *gs = 0x0;
622   if(!fGraphS){ 
623     fGraphS = new TObjArray(fNRefFigures);
624     fGraphS->SetOwner();
625     for(Int_t ig=0; ig<fNRefFigures; ig++){
626       gs = new TGraphErrors();
627       gs->SetLineColor(kRed);
628       gs->SetMarkerStyle(23);
629       gs->SetMarkerColor(kRed);
630       gs->SetNameTitle(Form("s_%d", ig), "");
631       fGraphS->AddAt(gs, ig);
632     }
633   }
634   if(!fGraphM){ 
635     fGraphM = new TObjArray(fNRefFigures);
636     fGraphM->SetOwner();
637     for(Int_t ig=0; ig<fNRefFigures; ig++){
638       gm = new TGraphErrors();
639       gm->SetLineColor(kBlue);
640       gm->SetMarkerStyle(7);
641       gm->SetMarkerColor(kBlue);
642       gm->SetNameTitle(Form("m_%d", ig), "");
643       fGraphM->AddAt(gm, ig);
644     }
645   }
646
647   TH2I *h2 = 0x0;
648   TH1D *h = 0x0;
649
650   // define models
651   TF1 f("f1", "gaus", -.5, .5);  
652
653   TF1 fb("fb", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]", -.5, .5);
654
655   TF1 fc("fc", "[0]*exp(-0.5*((x-[1])/[2])**2)+[3]*exp(-0.5*((x-[4])/[5])**2)", -.5, .5);
656
657   TCanvas *c = 0x0;
658   if(IsVisual()) c = new TCanvas("c", Form("%s Visual", GetName()), 500, 500);
659   char opt[5];
660   sprintf(opt, "%c%c", IsVerbose() ? ' ' : 'Q', IsVisual() ? ' ': 'N');
661
662
663   //PROCESS RESIDUAL DISTRIBUTIONS
664
665   // Clusters residuals
666   h2 = (TH2I *)(fContainer->At(kCluster));
667   gm = (TGraphErrors*)fGraphM->At(kCluster);
668   gs = (TGraphErrors*)fGraphS->At(kCluster);
669   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
670     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
671     h = h2->ProjectionY("py", ibin, ibin);
672     if(h->GetEntries()<100) continue;
673     AdjustF1(h, &f);
674
675     if(IsVisual()){c->cd(); c->SetLogy();}
676     h->Fit(&f, opt, "", -0.5, 0.5);
677     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
678     
679     Int_t ip = gm->GetN();
680     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
681     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
682     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
683     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
684   }
685
686   // Tracklet y residuals
687   h2 = (TH2I *)(fContainer->At(kTrackletY));
688   gm = (TGraphErrors*)fGraphM->At(kTrackletY);
689   gs = (TGraphErrors*)fGraphS->At(kTrackletY);
690   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
691     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
692     h = h2->ProjectionY("py", ibin, ibin);
693     if(h->GetEntries()<100) continue;
694     AdjustF1(h, &f);
695
696     if(IsVisual()){c->cd(); c->SetLogy();}
697     h->Fit(&f, opt, "", -0.5, 0.5);
698     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
699     
700     Int_t ip = gm->GetN();
701     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
702     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
703     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
704     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
705   }
706
707   // Tracklet phi residuals
708   h2 = (TH2I *)(fContainer->At(kTrackletPhi));
709   gm = (TGraphErrors*)fGraphM->At(kTrackletPhi);
710   gs = (TGraphErrors*)fGraphS->At(kTrackletPhi);
711   for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
712     Double_t phi = h2->GetXaxis()->GetBinCenter(ibin);
713     h = h2->ProjectionY("py", ibin, ibin);
714     if(h->GetEntries()<100) continue;
715     AdjustF1(h, &f);
716
717     if(IsVisual()){c->cd(); c->SetLogy();}
718     h->Fit(&f, opt, "", -0.5, 0.5);
719     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
720     
721     Int_t ip = gm->GetN();
722     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
723     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
724     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
725     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
726   }
727
728   if(!HasMCdata()){
729     if(c) delete c;
730     return kTRUE;
731   }
732
733   //PROCESS MC RESIDUAL DISTRIBUTIONS
734
735   // cluster y resolution
736   h2 = (TH2I*)fContainer->At(kMCcluster);
737   gm = (TGraphErrors*)fGraphM->At(kMCcluster);
738   gs = (TGraphErrors*)fGraphS->At(kMCcluster);
739   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
740     h = h2->ProjectionY("py", iphi, iphi);
741     if(h->GetEntries()<100) continue;
742     AdjustF1(h, &f);
743
744     if(IsVisual()){c->cd(); c->SetLogy();}
745     h->Fit(&f, opt, "", -0.5, 0.5);
746     if(IsVerbose()){
747       printf("phi[%d] mean[%e] sigma[%e]\n\n", iphi, 10.*f.GetParameter(1), 10.*f.GetParameter(2));
748     }
749     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
750
751     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
752     Int_t ip = gm->GetN();
753     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
754     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
755     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
756     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
757   }
758
759   // tracklet y resolution
760   h2 = (TH2I*)fContainer->At(kMCtrackletY);
761   gm = (TGraphErrors*)fGraphM->At(kMCtrackletY);
762   gs = (TGraphErrors*)fGraphS->At(kMCtrackletY);
763   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
764     h = h2->ProjectionY("py", iphi, iphi);
765     if(h->GetEntries()<100) continue;
766     AdjustF1(h, &f);
767
768     if(IsVisual()){c->cd(); c->SetLogy();}
769     h->Fit(&f, opt, "", -0.5, 0.5);
770     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
771
772     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
773     Int_t ip = gm->GetN();
774     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
775     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
776     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
777     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
778   }
779
780   // tracklet z resolution
781   h2 = (TH2I*)fContainer->At(kMCtrackletZ);
782   gm = (TGraphErrors*)fGraphM->At(kMCtrackletZ);
783   gs = (TGraphErrors*)fGraphS->At(kMCtrackletZ);
784   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
785     h = h2->ProjectionY("py", iphi, iphi);
786     if(h->GetEntries()<100) continue;
787     AdjustF1(h, &fb);
788
789     if(IsVisual()){c->cd(); c->SetLogy();}
790     h->Fit(&fb, opt, "", -0.5, 0.5);
791     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
792
793     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
794     Int_t ip = gm->GetN();
795     gm->SetPoint(ip, phi, 10.*fb.GetParameter(1));
796     gm->SetPointError(ip, 0., 10.*fb.GetParError(1));
797     gs->SetPoint(ip, phi, 10.*fb.GetParameter(2));
798     gs->SetPointError(ip, 0., 10.*fb.GetParError(2));
799   }
800
801   //tracklet phi resolution
802   h2 = (TH2I*)fContainer->At(kMCtrackletPhi);
803   gm = (TGraphErrors*)fGraphM->At(kMCtrackletPhi);
804   gs = (TGraphErrors*)fGraphS->At(kMCtrackletPhi);
805   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
806     h = h2->ProjectionY("py", iphi, iphi);
807     if(h->GetEntries()<100) continue;
808
809     if(IsVisual()){c->cd(); c->SetLogy();}
810     h->Fit(&f, opt, "", -0.5, 0.5);
811     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
812
813     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
814     Int_t ip = gm->GetN();
815     gm->SetPoint(ip, phi, f.GetParameter(1));
816     gm->SetPointError(ip, 0., f.GetParError(1));
817     gs->SetPoint(ip, phi, f.GetParameter(2));
818     gs->SetPointError(ip, 0., f.GetParError(2));
819   }
820
821   // track y resolution
822   h2 = (TH2I*)fContainer->At(kMCtrackY);
823   gm = (TGraphErrors*)fGraphM->At(kMCtrackY);
824   gs = (TGraphErrors*)fGraphS->At(kMCtrackY);
825   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
826     h = h2->ProjectionY("py", iphi, iphi);
827     if(h->GetEntries()<100) continue;
828     AdjustF1(h, &f);
829
830     if(IsVisual()){c->cd(); c->SetLogy();}
831     h->Fit(&f, opt, "", -0.5, 0.5);
832     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
833
834     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
835     Int_t ip = gm->GetN();
836     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
837     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
838     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
839     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
840   }
841
842   // track z resolution
843   h2 = (TH2I*)fContainer->At(kMCtrackZ);
844   gm = (TGraphErrors*)fGraphM->At(kMCtrackZ);
845   gs = (TGraphErrors*)fGraphS->At(kMCtrackZ);
846   for(Int_t iphi=1; iphi<=h2->GetNbinsX(); iphi++){
847     h = h2->ProjectionY("pz", iphi, iphi);
848     if(h->GetEntries()<70) continue;
849     AdjustF1(h, &f);
850
851     if(IsVisual()){c->cd(); c->SetLogy();}
852     h->Fit(&f, opt, "", -0.5, 0.5);
853     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
854
855     Double_t phi = h2->GetXaxis()->GetBinCenter(iphi);
856     Int_t ip = gm->GetN();
857     gm->SetPoint(ip, phi, 10.*f.GetParameter(1));
858     gm->SetPointError(ip, 0., 10.*f.GetParError(1));
859     gs->SetPoint(ip, phi, 10.*f.GetParameter(2));
860     gs->SetPointError(ip, 0., 10.*f.GetParError(2));
861   }
862
863   // track Pt resolution
864   h2 = (TH2I*)fContainer->At(kMCtrackPt);
865   TAxis *ax = h2->GetXaxis();
866   gm = (TGraphErrors*)fGraphM->At(kMCtrackPt);
867   gs = (TGraphErrors*)fGraphS->At(kMCtrackPt);
868   TF1 fg("fg", "gaus", -1.5, 1.5);
869   TF1 fl("fl", "landau", -4., 15.);
870   TF1 fgl("fgl", "gaus(0)+landau(3)", -5., 20.);
871   for(Int_t ip=1; ip<=ax->GetNbins(); ip++){
872     h = h2->ProjectionY("ppt", ip, ip);
873     if(h->GetEntries()<70) continue;
874
875     h->Fit(&fg, "QN", "", -1.5, 1.5);
876     fgl.SetParameter(0, fg.GetParameter(0));
877     fgl.SetParameter(1, fg.GetParameter(1));
878     fgl.SetParameter(2, fg.GetParameter(2));
879     h->Fit(&fl, "QN", "", -4., 15.);
880     fgl.SetParameter(3, fl.GetParameter(0));
881     fgl.SetParameter(4, fl.GetParameter(1));
882     fgl.SetParameter(5, fl.GetParameter(2));
883
884     if(IsVisual()){c->cd(); c->SetLogy();}
885     h->Fit(&fgl, opt, "", -5., 20.);
886     if(IsVisual()){c->Modified(); c->Update(); gSystem->Sleep(500);}
887
888     Float_t invpt = ax->GetBinCenter(ip);
889     Int_t ip = gm->GetN();
890     gm->SetPoint(ip, invpt, fgl.GetParameter(1));
891     gm->SetPointError(ip, 0., fgl.GetParError(1));
892     gs->SetPoint(ip, invpt, fgl.GetParameter(2)*invpt);
893     gs->SetPointError(ip, 0., fgl.GetParError(2));
894     // fgl.GetParameter(4) // Landau MPV
895     // fgl.GetParameter(5) // Landau Sigma
896   }
897
898
899   if(c) delete c;
900   return kTRUE;
901 }
902
903
904 //________________________________________________________
905 void AliTRDtrackingResolution::Terminate(Option_t *)
906 {
907   if(fDebugStream){ 
908     delete fDebugStream;
909     fDebugStream = 0x0;
910     fDebugLevel = 0;
911   }
912   if(HasPostProcess()) PostProcess();
913 }
914
915 //________________________________________________________
916 void AliTRDtrackingResolution::AdjustF1(TH1 *h, TF1 *f)
917 {
918 // Helper function to avoid duplication of code
919 // Make first guesses on the fit parameters
920
921   // find the intial parameters of the fit !! (thanks George)
922   Int_t nbinsy = Int_t(.5*h->GetNbinsX());
923   Double_t sum = 0.;
924   for(Int_t jbin=nbinsy-4; jbin<=nbinsy+4; jbin++) sum+=h->GetBinContent(jbin); sum/=9.;
925   f->SetParLimits(0, 0., 3.*sum);
926   f->SetParameter(0, .9*sum);
927
928   f->SetParLimits(1, -.2, .2);
929   f->SetParameter(1, -0.1);
930
931   f->SetParLimits(2, 0., 4.e-1);
932   f->SetParameter(2, 2.e-2);
933   if(f->GetNpar() <= 4) return;
934
935   f->SetParLimits(3, 0., sum);
936   f->SetParameter(3, .1*sum);
937
938   f->SetParLimits(4, -.3, .3);
939   f->SetParameter(4, 0.);
940
941   f->SetParLimits(5, 0., 1.e2);
942   f->SetParameter(5, 2.e-1);
943 }
944
945 //________________________________________________________
946 TObjArray* AliTRDtrackingResolution::Histos()
947 {
948   if(fContainer) return fContainer;
949
950   fContainer  = new TObjArray(10);
951   //fContainer->SetOwner(kTRUE);
952
953   TH1 *h = 0x0;
954   // cluster to tracklet residuals [2]
955   if(!(h = (TH2I*)gROOT->FindObject("hCl"))){
956     h = new TH2I("hCl", "Clusters-Track Residuals", 21, -.33, .33, 100, -.5, .5);
957     h->GetXaxis()->SetTitle("tg(#phi)");
958     h->GetYaxis()->SetTitle("#Delta y [cm]");
959     h->GetZaxis()->SetTitle("entries");
960   } else h->Reset();
961   fContainer->AddAt(h, kCluster);
962
963   // tracklet to track residuals [2]
964   if(!(h = (TH2I*)gROOT->FindObject("hTrkltY"))){
965     h = new TH2I("hTrkltY", "Tracklets-Track Residuals (Y)", 21, -.33, .33, 100, -.5, .5);
966     h->GetXaxis()->SetTitle("tg(#phi)");
967     h->GetYaxis()->SetTitle("#Delta y [cm]");
968     h->GetZaxis()->SetTitle("entries");
969   } else h->Reset();
970   fContainer->AddAt(h, kTrackletY);
971
972   // tracklet to track residuals angular [2]
973   if(!(h = (TH2I*)gROOT->FindObject("hTrkltPhi"))){
974     h = new TH2I("hTrkltPhi", "Tracklets-Track Residuals (#Phi)", 21, -.33, .33, 100, -.5, .5);
975     h->GetXaxis()->SetTitle("tg(#phi)");
976     h->GetYaxis()->SetTitle("#Delta phi [#circ]");
977     h->GetZaxis()->SetTitle("entries");
978   } else h->Reset();
979   fContainer->AddAt(h, kTrackletPhi);
980
981
982   // Resolution histos
983   if(!HasMCdata()) return fContainer;
984
985   // cluster y resolution [0]
986   if(!(h = (TH2I*)gROOT->FindObject("hMCcl"))){
987     h = new TH2I("hMCcl", "Cluster Resolution", 31, -.48, .48, 100, -.5, .5);
988     h->GetXaxis()->SetTitle("tg(#phi)");
989     h->GetYaxis()->SetTitle("#Delta y [cm]");
990     h->GetZaxis()->SetTitle("entries");
991   } else h->Reset();
992   fContainer->AddAt(h, kMCcluster);
993
994   // tracklet y resolution [0]
995   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkltY"))){
996     h = new TH2I("hMCtrkltY", "Tracklet Resolution (Y)", 31, -.48, .48, 100, -.5, .5);
997     h->GetXaxis()->SetTitle("tg(#phi)");
998     h->GetYaxis()->SetTitle("#Delta y [cm]");
999     h->GetZaxis()->SetTitle("entries");
1000   } else h->Reset();
1001   fContainer->AddAt(h, kMCtrackletY);
1002
1003   // tracklet y resolution [0]
1004   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkltZ"))){
1005     h = new TH2I("hMCtrkltZ", "Tracklet Resolution (Z)", 31, -.48, .48, 100, -.5, .5);
1006     h->GetXaxis()->SetTitle("tg(#theta)");
1007     h->GetYaxis()->SetTitle("#Delta z [cm]");
1008     h->GetZaxis()->SetTitle("entries");
1009   } else h->Reset();
1010   fContainer->AddAt(h, kMCtrackletZ);
1011
1012   // tracklet angular resolution [1]
1013   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkltPhi"))){
1014     h = new TH2I("hMCtrkltPhi", "Tracklet Resolution (#Phi)", 31, -.48, .48, 100, -10., 10.);
1015     h->GetXaxis()->SetTitle("tg(#phi)");
1016     h->GetYaxis()->SetTitle("#Delta #phi [deg]");
1017     h->GetZaxis()->SetTitle("entries");
1018   } else h->Reset();
1019   fContainer->AddAt(h, kMCtrackletPhi);
1020
1021   // Kalman track y resolution
1022   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkY"))){
1023     h = new TH2I("hMCtrkY", "Kalman Track Resolution (Y)", 31, -.48, .48, 100, -.5, .5);
1024     h->GetXaxis()->SetTitle("tg(#phi)");
1025     h->GetYaxis()->SetTitle("#Delta y [cm]");
1026     h->GetZaxis()->SetTitle("entries");
1027   } else h->Reset();
1028   fContainer->AddAt(h, kMCtrackY);
1029
1030   // Kalman track Z resolution
1031   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkZ"))){
1032     h = new TH2I("hMCtrkZ", "Kalman Track Resolution (Z)", 20, -1., 1., 100, -1.5, 1.5);
1033     h->GetXaxis()->SetTitle("tg(#theta)");
1034     h->GetYaxis()->SetTitle("#Delta z [cm]");
1035     h->GetZaxis()->SetTitle("entries");
1036   } else h->Reset();
1037   fContainer->AddAt(h, kMCtrackZ);
1038
1039   // Kalman track Pt resolution
1040   if(!(h = (TH2I*)gROOT->FindObject("hMCtrkPt"))){
1041     h = new TH2I("hMCtrkPt", "Kalman Track Resolution (Pt)", 100, 0., 2., 150, -5., 20.);
1042     h->GetXaxis()->SetTitle("1/p_{t}");
1043     h->GetYaxis()->SetTitle("#Delta p_{t} [GeV/c]");
1044     h->GetZaxis()->SetTitle("entries");
1045   } else h->Reset();
1046   fContainer->AddAt(h, kMCtrackPt);
1047
1048   return fContainer;
1049 }
1050
1051
1052 //________________________________________________________
1053 void AliTRDtrackingResolution::SetRecoParam(AliTRDrecoParam *r)
1054 {
1055
1056   fReconstructor->SetRecoParam(r);
1057 }