ROOT logo
/**************************************************************************
* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
*                                                                        *
* Author: The ALICE Off-line Project.                                    *
* Contributors are mentioned in the code where appropriate.              *
*                                                                        *
* Permission to use, copy, modify and distribute this software and its   *
* documentation strictly for non-commercialf purposes is hereby granted   *
* without fee, provided that the above copyright notice appears in all   *
* copies and that both the copyright notice and this permission notice   *
* appear in the supporting documentation. The authors make no claims     *
* about the suitability of this software for any purpose. It is          *
* provided "as is" without express or implied warranty.                  *
**************************************************************************/

/* $Id: AliTRDclusterResolution.cxx */

///////////////////////////////////////////////////////////////////////////////
//                                                                           //
//  TRD cluster error parameterization                                        //
//                                                                           //
// This class is designed to produce the reference plots for a detailed study//
// and parameterization of TRD cluster errors. The following effects are taken//
// into account :                                                            //
//   - dependence with the total charge of the cluster                       //
//   - dependence with the distance from the center pad. This is monitored 
// for each layer individually since the pad size varies with layer
//   - dependence with the drift length - here the influence of anisochronity 
// and diffusion are searched
//   - dependence with the distance to the anode wire - anisochronity effects
//   - dependence with track angle (for y resolution)
// The correlation between effects is taken into account. 
// 
// Since magnetic field plays a very important role in the TRD measurement 
// the ExB correction is forced by the setter function SetExB(Int_t). The 
// argument is the detector index, if none is specified all will be 
// considered.
// 
// Two cases are of big importance.
//   - comparison with MC
//   - comparison with Kalman fit. In this case the covariance matrix of the
// Kalman fit are needed.
// 
// The functionalities implemented in this class are based on the storage 
// class AliTRDclusterInfo.
// 
// The Method
// ----------
// 
// The method to disentangle s_y and s_x is based on the relation (see also fig.)
// BEGIN_LATEX
// #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}})
// END_LATEX
// with
// BEGIN_LATEX
// #sigma^{2}_{x_{c}} #approx 0 
// END_LATEX
// we suppose the chamber is well calibrated for t_{0} and aligned in
// radial direction. 
//
// Clusters can be radially shifted due to three causes:
//   - globally shifted - due to residual misalignment/miscalibration(t0)
//   - locally shifted - due to different local drift velocity from the mean
//   - randomly shifted - due to neighboring (radial direction) clusters 
// charge induced by asymmetry of the TRF.
//
// We estimate this effects by the relations:
// BEGIN_LATEX
// #mu_{y} = tg(#alpha_{L})*#Delta x_{d}(...) + tg(#phi-#alpha_{L})*(#Delta x_{c}(...) + #Delta x_{d}(...))
// END_LATEX
// where
// BEGIN_LATEX
// #Delta x_{d}(...) = (<v_{d}> + #delta v_{d}(x_{d}, d)) * (t + t^{*}(Q))
// END_LATEX
// and we specified explicitely the variation of drift velocity parallel 
// with the track (x_{d}) and perpendicular to it due to anisochronity (d).
// 
// For estimating the contribution from asymmetry of TRF the following
// parameterization is being used
// BEGIN_LATEX
// t^{*}(Q) = #delta_{0} * #frac{Q_{t+1} - Q_{t-1}}{Q_{t-1} + Q_{t} + Q_{t+1}}
// END_LATEX
//
//
// Clusters can also be r-phi shifted due to:
//   - wrong PRF or wrong cuts at digits level
//The following correction is applied :
// BEGIN_LATEX
// <#Delta y> = a + b * sin(c*y_{pw})
// END_LATEX

// The Models
//
//   Parameterization against total charge
//
// Obtained for B=0T at phi=0. All other effects integrated out.
// BEGIN_LATEX
// #sigma^{2}_{y}(Q) = #sigma^{2}_{y}(...) + b(#frac{1}{Q} - #frac{1}{Q_{0}}) 
// END_LATEX
// For B diff 0T the error of the average ExB correction error has to be subtracted !! 
//
//   Parameterization Sx
//
// The parameterization of the error in the x direction can be written as
// BEGIN_LATEX
// #sigma_{x} = #sigma_{x}^{||} + #sigma_{x}^{#perp}
// END_LATEX
//
// where the parallel component is given mainly by the TRF width while 
// the perpendicular component by the anisochronity. The model employed for 
// the parallel is gaus(0)+expo(3) with the following parameters
// 1  C   5.49018e-01   1.23854e+00   3.84540e-04  -8.21084e-06
// 2  M   7.82999e-01   6.22531e-01   2.71272e-04  -6.88485e-05
// 3  S   2.74451e-01   1.13815e+00   2.90667e-04   1.13493e-05
// 4  E1  2.53596e-01   1.08646e+00   9.95591e-05  -2.11625e-05
// 5  E2 -2.40078e-02   4.26520e-01   4.67153e-05  -2.35392e-04
//
// and perpendicular to the track is pol2 with the parameters
//
// Par_0 = 0.190676 +/- 0.41785
// Par_1 = -3.9269  +/- 7.49862
// Par_2 = 14.7851  +/- 27.8012
//
//   Parameterization Sy
//
// The parameterization of the error in the y direction along track uses
// BEGIN_LATEX
// #sigma_{y}^{||} = #sigma_{y}^{0} -a*exp(1/(x-b))
// END_LATEX
//
// with following values for the parameters:
// 1  sy0 2.60967e-01   2.99652e-03   7.82902e-06  -1.89636e-04
// 2  a  -7.68941e+00   1.87883e+00   3.84539e-04   9.38268e-07
// 3  b  -3.41160e-01   7.72850e-02   1.63231e-05   2.51602e-05
//
//==========================================================================
// Example how to retrive reference plots from the task
// void steerClErrParam(Int_t fig=0)
// {
//   gSystem->Load("libANALYSIS.so");
//   gSystem->Load("libTRDqaRec.so");
// 
//   // initialize DB manager
//   AliCDBManager *cdb = AliCDBManager::Instance();
//   cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
//   cdb->SetRun(0);
//   // initialize magnetic field.
//   AliMagFCheb *field=new AliMagFCheb("Maps","Maps", 2, 1., 10., AliMagFCheb::k5kG);
//   AliTracker::SetFieldMap(field, kTRUE);
// 
//   AliTRDclusterResolution *res = new AliTRDclusterResolution();
//   res->SetMCdata();
//   res->Load("TRD.TaskClErrParam.root");
//   res->SetExB();  
//   res->SetVisual(); 
//   //res->SetSaveAs();
//   res->SetProcessCharge(kFALSE);
//   res->SetProcessCenterPad(kFALSE);
//   //res->SetProcessMean(kFALSE);
//   res->SetProcessSigma(kFALSE);
//   if(!res->PostProcess()) return;
//   new TCanvas;
//   res->GetRefFigure(fig);
// }
//
//  Authors:                                                              //
//    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
////////////////////////////////////////////////////////////////////////////

#include "AliTRDclusterResolution.h"
#include "AliTRDresolution.h"
#include "AliTRDinfoGen.h"
#include "info/AliTRDclusterInfo.h"
#include "info/AliTRDeventInfo.h"

#include "AliTRDcalibDB.h"
#include "Cal/AliTRDCalROC.h"
#include "Cal/AliTRDCalDet.h"
#include "AliTRDCommonParam.h"
#include "AliTRDgeometry.h"
#include "AliTRDpadPlane.h"
#include "AliTRDcluster.h"
#include "AliTRDseedV1.h"

#include "AliESDEvent.h"
#include "AliCDBManager.h"

#include "TROOT.h"
#include "TSystem.h"
#include "TMath.h"
#include "TLinearFitter.h"
#include "TGeoGlobalMagField.h"
#include <TGeoMatrix.h>
#include "TObjArray.h"
#include "TTree.h"
#include "TH2I.h"
#include "TH3S.h"
#include "TAxis.h"
#include "TF1.h"
#include "TCanvas.h"
#include "TLegend.h"
#include "TGraphErrors.h"
#include "TLine.h"

ClassImp(AliTRDclusterResolution)

const Float_t AliTRDclusterResolution::fgkTimeBinLength = 1./ AliTRDCommonParam::Instance()->GetSamplingFrequency();
//_______________________________________________________
AliTRDclusterResolution::AliTRDclusterResolution()
  : AliTRDrecoTask()
  ,fCanvas(NULL)
  ,fInfo(NULL)
  ,fResults(NULL)
  ,fSubTaskMap(0)
  ,fUseCalib(7)
  ,fDet(-1)
  ,fCol(-1)
  ,fRow(-1)
  ,fExB(0.)
  ,fDt(0.)
  ,fDl(0.)
  ,fVdrift(1.5)
  ,fT0(0.)
  ,fGain(1.)
  ,fXch(0.)
  ,fZch(0.)
  ,fH(0.)
  ,fDyRange(1.3)
  ,fLy(0)
  ,fT(0.)
  ,fX(0.)
  ,fY(0.)
  ,fZ(0.)
{
// Constructor
  SetNameTitle("ClErrCalib", "Cluster Error Parameterization");
  memset(fR, 0, 4*sizeof(Float_t));
  memset(fP, 0, 4*sizeof(Float_t));
}

//_______________________________________________________
AliTRDclusterResolution::AliTRDclusterResolution(const char *name)
  : AliTRDrecoTask(name, "Cluster Error Parameterization")
  ,fCanvas(NULL)
  ,fInfo(NULL)
  ,fResults(NULL)
  ,fSubTaskMap(0)
  ,fUseCalib(7)
  ,fDet(-1)
  ,fCol(-1)
  ,fRow(-1)
  ,fExB(0.)
  ,fDt(0.)
  ,fDl(0.)
  ,fVdrift(1.5)
  ,fT0(0.)
  ,fGain(1.)
  ,fXch(0.)
  ,fZch(0.)
  ,fH(0.)
  ,fDyRange(1.3)
  ,fLy(0)
  ,fT(0.)
  ,fX(0.)
  ,fY(0.)
  ,fZ(0.)
{
// Constructor

  memset(fR, 0, 4*sizeof(Float_t));
  memset(fP, 0, 4*sizeof(Float_t));

  // By default register all analysis
  // The user can switch them off in his steering macro
  SetProcess(kYRes);
  SetProcess(kYSys);
  SetProcess(kMean);
  SetProcess(kSigm);
}

//_______________________________________________________
AliTRDclusterResolution::~AliTRDclusterResolution()
{
// Destructor

  if(fCanvas) delete fCanvas;
  if(fResults){
    fResults->Delete();
    delete fResults;
  }
}

//_______________________________________________________
void AliTRDclusterResolution::UserCreateOutputObjects()
{
// Build and post histo container.
// Actual population of the container with histo is done in function Histos.

  if(!fContainer) fContainer = new TObjArray(kNtasks);
 //fContainer->SetOwner(kTRUE);
  PostData(1, fContainer);
}

//_______________________________________________________
Bool_t AliTRDclusterResolution::GetRefFigure(Int_t ifig)
{
// Steering function to retrieve performance plots

  if(!fResults) return kFALSE;
  TLegend *leg = NULL;
  TList *l = NULL;
  TObjArray *arr = NULL;
  TTree *t = NULL;
  TH2 *h2 = NULL;TH1 *h1 = NULL;
  TGraphErrors *gm(NULL), *gs(NULL), *gp(NULL);
  switch(ifig){
  case kYRes:
    if(!(arr = (TObjArray*)fResults->At(kYRes))) break;
    if(!(gm = (TGraphErrors*)arr->At(0))) break;
    if(!(gs = (TGraphErrors*)arr->At(1))) break;
    if(!(gp = (TGraphErrors*)arr->At(2))) break;
    leg= new TLegend(.7, .7, .9, .95);
    leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetFillStyle(0);
    gs->Draw("apl"); leg->AddEntry(gs, "Sigma / Resolution", "pl");
    gs->GetHistogram()->GetYaxis()->SetRangeUser(-50., 700.);
    gs->GetHistogram()->SetXTitle("Q [a.u.]");
    gs->GetHistogram()->SetYTitle("y - x tg(#alpha_{L}) [#mum]");
    gm->Draw("pl");leg->AddEntry(gm, "Mean / Systematics", "pl");
    gp->Draw("pl");leg->AddEntry(gp, "Abundance / Probability", "pl");
    leg->Draw();
    return kTRUE;
  case kYSys:
    if(!(arr = (TObjArray*)fResults->At(kYSys))) break;
    gPad->Divide(2, 1); l = gPad->GetListOfPrimitives();
    ((TVirtualPad*)l->At(0))->cd();
    ((TTree*)arr->At(0))->Draw(Form("y:t>>h(%d, -0.5, %f, 51, -.51, .51)",AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5),
            "m[0]*(ly==0&&abs(m[0])<1.e-1)", "colz");
    ((TVirtualPad*)l->At(1))->cd();
    leg= new TLegend(.7, .7, .9, .95);
    leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetFillStyle(0);
    leg->SetHeader("TRD Plane"); 
    for(Int_t il = 1; il<=AliTRDgeometry::kNlayer; il++){
      if(!(gm = (TGraphErrors*)arr->At(il))) return kFALSE;
      gm->Draw(il>1?"pc":"apc"); leg->AddEntry(gm, Form("%d", il-1), "pl");
      if(il>1) continue;
      gm->GetHistogram()->SetXTitle("t_{drift} [tb]");
      gm->GetHistogram()->SetYTitle("#sigma_{y}(x|cen=0) [#mum]");
      gm->GetHistogram()->GetYaxis()->SetRangeUser(150., 500.);
    }
    leg->Draw();
    return kTRUE;
  case kSigm:
    if(!(t = (TTree*)fResults->At(kSigm))) break;
    t->Draw("z:t>>h2x(23, 0.1, 2.4, 25, 0., 2.5)","sx*(1)", "lego2fb");
    h2 = (TH2F*)gROOT->FindObject("h2x");
    printf("  const Double_t sx[24][25]={\n");
    for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
      printf("    {");
      for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
        printf("%6.4f ", h2->GetBinContent(ix, iy));
      }
      printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
    }
    printf("  };\n");
    gPad->Divide(2, 1, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
    ((TVirtualPad*)l->At(0))->cd();
    h1 = h2->ProjectionX("hsx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
    h1->SetYTitle("<#sigma_{x}> [#mum]");
    h1->SetXTitle("t_{drift} [#mus]");
    h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); h1->Draw("pc");

    t->Draw("z:t>>h2y(23, 0.1, 2.4, 25, 0., 2.5)","sy*(1)", "lego2fb");
    h2 = (TH2F*)gROOT->FindObject("h2y");
    printf("  const Double_t sy[24][25]={\n");
    for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
      printf("    {");
      for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
        printf("%6.4f ", h2->GetBinContent(ix, iy));
      }
      printf("%6.4f},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
    }
    printf("  };\n");
    ((TVirtualPad*)l->At(1))->cd();
    h1 = h2->ProjectionX("hsy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
    h1->SetYTitle("<#sigma_{y}> [#mum]");
    h1->SetXTitle("t_{drift} [#mus]");
    h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); h1->Draw("pc");
    return kTRUE;
  case kMean:
    if(!(t = (TTree*)fResults->At(kMean))) break;
    if(!t->Draw(Form("z:t>>h2x(%d, -0.5, %3.1f, %d, 0., 2.5)", 
      AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5, kND),
      "dx*(1)", "goff")) break;
    h2 = (TH2F*)gROOT->FindObject("h2x");
    printf("  const Double_t dx[%d][%d]={\n", AliTRDseedV1::kNtb, kND);
    for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
      printf("    {");
      for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
        printf("%+6.4e, ", h2->GetBinContent(ix, iy));
      }
      printf("%+6.4e},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
    }
    printf("  };\n");
    gPad->Divide(2, 2, 1.e-5, 1.e-5); l = gPad->GetListOfPrimitives();
    ((TVirtualPad*)l->At(0))->cd();
    h2->Draw("lego2fb");
    ((TVirtualPad*)l->At(2))->cd();
    h1 = h2->ProjectionX("hdx_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
    h1->SetYTitle("<#deltax> [#mum]");
    h1->SetXTitle("t_{drift} [tb]");
    //h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); 
    h1->Draw("pc");

    if(!t->Draw(Form("z:t>>h2y(%d, -0.5, %3.1f, %d, 0., 2.5)", 
      AliTRDseedV1::kNtb, AliTRDseedV1::kNtb-0.5, kND),
      "dy*(1)", "goff")) break;
    h2 = (TH2F*)gROOT->FindObject("h2y");
    printf("  const Double_t dy[%d][%d]={\n", AliTRDseedV1::kNtb, kND);
    for(Int_t ix=1; ix<=h2->GetNbinsX(); ix++){
      printf("    {");
      for(Int_t iy=1; iy<h2->GetNbinsY(); iy++){
        printf("%+6.4e ", h2->GetBinContent(ix, iy));
      }
      printf("%+6.4e},\n", h2->GetBinContent(ix, h2->GetNbinsY()));
    }
    printf("  };\n");
    ((TVirtualPad*)l->At(1))->cd();
    h2->Draw("lego2fb");
    ((TVirtualPad*)l->At(3))->cd();
    h1 = h2->ProjectionX("hdy_pxx"); h1->Scale(1.e4/kND); h1->SetMarkerStyle(24);
    h1->SetYTitle("<#deltay> [#mum]");
    h1->SetXTitle("t_{drift} [tb]");
    //h1->GetXaxis()->SetRange(2, AliTRDseedV1::kNtb-1); 
    h1->Draw("pc");

    return kTRUE;
  default:
    break;
  }
  AliWarning("No container/data found.");
  return kFALSE;
}

//_______________________________________________________
TObjArray* AliTRDclusterResolution::Histos()
{
// Retrieve histograms array if already build or build it

  if(!fContainer){
    fContainer = new TObjArray(kNtasks);
    //fContainer->SetOwner(kTRUE);
  }
  if(fContainer->GetEntries() == kNtasks) return fContainer;

  TH3S *h3(NULL);TH2I *h2(NULL);
  TObjArray *arr(NULL);
  if(!HasGlobalPosition() && !LoadGlobalChamberPosition()) return NULL;
  Float_t tgt(fZch/fXch), htgt(fH*tgt);
  
  // SYSTEMATIC PLOTS
  fContainer->AddAt(arr = new TObjArray(3), kYSys);
  arr->SetName("SysY");
  // systematic plot on pw and q (dydx=ExB+h*dzdx)
  if(!(h3=(TH3S*)gROOT->FindObject(Form("Sys%s%03d", (HasMCdata()?"MC":"") ,fDet)))) {
    h3 = new TH3S(
      Form("Sys%s%03d", (HasMCdata()?"MC":""),fDet),
      Form(" Det[%d] Col[%d] Row[%d];log q [a.u.];#deltay [pw];#Delta y[cm]", fDet, fCol, fRow),
      45,   2., 6.5,            // log(q) [a.u.]
      25, -.51, .51,            // y [pw]
      60, -fDyRange, fDyRange); // dy [cm]
  } h3->Reset();
  arr->AddAt(h3, 0);
  // systematic plot on tb (only for dydx = h*tgt + exb and MPV q)
  if(!(h2 = (TH2I*)gROOT->FindObject(Form("SysTb%s%03d", (HasMCdata()?"MC":""), fDet)))){
    h2 = new TH2I(Form("SysTb%s%03d", (HasMCdata()?"MC":""), fDet),
    Form(" Det[%d] Col[%d] Row[%d];t [time bin];#Delta y[cm]", fDet, fCol, fRow),
    AliTRDseedV1::kNtb, -.5, AliTRDseedV1::kNtb-0.5,   // t [tb]
    60, -fDyRange, fDyRange);                          // dy [cm]
  } h2->Reset();
  arr->AddAt(h2, 1);
  // systematic plot on tgp and tb (for MPV q)
  if(!(h3=(TH3S*)gROOT->FindObject(Form("SysTbTgp%s%03d", (HasMCdata()?"MC":""), fDet)))){
    h3 = new TH3S(
      Form("SysTbTgp%s%03d", (HasMCdata()?"MC":""), fDet),
      Form(" Det[%d];t [time bin];tg(#phi) - h*tg(#theta) %s;#Delta y[cm]", fDet, fExB>1.e-5?"- tg(#alpha_{L})":""),
      AliTRDseedV1::kNtb, -.5, AliTRDseedV1::kNtb-0.5,   // t [tb]
      36, fExB-.18, fExB+.18,                            // tgp-h tgt-tg(aL)
      60, -fDyRange, fDyRange);                          // dy
  } h3->Reset();
  arr->AddAt(h3, 2);

  //  RESOLUTION/PULLS PLOTS
  fContainer->AddAt(arr = new TObjArray(6), kYRes);
  arr->SetName("ResY");
  // resolution plot on pw and q (for dydx=0 && B=0) np = 3 and for tb in [5, 20]
  TObjArray *arrt(NULL);
  arr->AddAt(arrt = new TObjArray(16), 0);
  arrt->SetName("PwQvsX");
  for(Int_t it(0); it<=15; it++){
    if(!(h3=(TH3S*)gROOT->FindObject(Form("Res%s%03d%02d", (HasMCdata()?"MC":"") ,fDet, it)))) {
      h3 = new TH3S(
        Form("Res%s%03d%02d", (HasMCdata()?"MC":""),fDet, it),
        Form(" Det[%d] TB[%d];log q [a.u];#deltay [pw];#Delta y[cm]", fDet, it+5),
        4,   2., 6.,            // log(q) [a.u]
        5, -.51, .51,            // y [pw]
        60, -fDyRange, fDyRange); // dy
    } h3->Reset();
    arrt->AddAt(h3, it);
  }
  // Pull plot on pw and q (for dydx=0 && B=0)
  if(!(h3=(TH3S*)gROOT->FindObject(Form("Pull%s%03d", (HasMCdata()?"MC":""), fDet)))){
    h3 = new TH3S(
      Form("Pull%s%03d", (HasMCdata()?"MC":""), fDet),
      Form(" Det[%d] Col[%d] Row[%d];log q [a.u.];#deltay [pw];#Delta y[cm]/#sigma_{y}", fDet, fCol, fRow),
      4,   2., 6.,            // log(q) [a.u]
      5, -.51, .51,            // y [pw]
      60, -4., 4.);             // dy/sy
  } h3->Reset();
  arr->AddAt(h3, 1);
  // resolution/pull plot on tb (for dydx=0 && B=0 && MPV q)
  if(!(h3 = (TH3S*)gROOT->FindObject(Form("ResPullTb%s%03d", (HasMCdata()?"MC":""), fDet)))){
    h3 = new TH3S(Form("ResPullTb%s%03d", (HasMCdata()?"MC":""), fDet),
    Form(" Det[%d] Col[%d] Row[%d];t [time bin];#Delta y[cm];#Delta y/#sigma_{y}", fDet, fCol, fRow),
    AliTRDseedV1::kNtb, -.5, AliTRDseedV1::kNtb-0.5,   // t [tb]
    60, -fDyRange, fDyRange,                           // dy [cm]
    60, -4., 4.);                                      // dy/sy
  } h3->Reset();
  arr->AddAt(h3, 2);
  // resolution plot on pw and q (for dydx=0 && B=0) np = 2
  if(!(h3=(TH3S*)gROOT->FindObject(Form("Res2%s%03d", (HasMCdata()?"MC":"") ,fDet)))) {
    h3 = new TH3S(
      Form("Res2%s%03d", (HasMCdata()?"MC":""),fDet),
      Form(" Det[%d] Col[%d] Row[%d];log q [a.u];#deltay [pw];#Delta y[cm]", fDet, fCol, fRow),
      4,   2., 6.,            // log(q) [a.u]
      5, -.51, .51,            // y [pw]
      60, -fDyRange, fDyRange); // dy
  } h3->Reset();
  arr->AddAt(h3, 3);
  // resolution plot on pw and q (for dydx=0 && B=0) np = 4
  if(!(h3=(TH3S*)gROOT->FindObject(Form("Res4%s%03d", (HasMCdata()?"MC":"") ,fDet)))) {
    h3 = new TH3S(
      Form("Res4%s%03d", (HasMCdata()?"MC":""),fDet),
      Form(" Det[%d] Col[%d] Row[%d];log q [a.u];#deltay [pw];#Delta y[cm]", fDet, fCol, fRow),
      4,   2., 6.,            // log(q) [a.u]
      5, -.51, .51,            // y [pw]
      60, -fDyRange, fDyRange); // dy
  } h3->Reset();
  arr->AddAt(h3, 4);
  // systemtic plot of tb on pw and q (for dydx=0 && B=0)
  if(!(h3=(TH3S*)gROOT->FindObject(Form("SysTbPwQ%s%03d", (HasMCdata()?"MC":"") ,fDet)))) {
    h3 = new TH3S(
      Form("SysTbPwQ%s%03d", (HasMCdata()?"MC":""),fDet),
      Form(" Det[%d] Col[%d] Row[%d];log q [a.u];#deltay [pw];t [time bin]", fDet, fCol, fRow),
      4,   2., 6.,            // log(q) [a.u]
      5, -.51, .51,            // y [pw]
      AliTRDseedV1::kNtb, -.5, AliTRDseedV1::kNtb-0.5);   // t [tb]
  } h3->Reset();
  arr->AddAt(h3, 5);



  fContainer->AddAt(arr = new TObjArray(AliTRDseedV1::kNtb), kSigm);
  arr->SetName("Resolution");
  for(Int_t it=0; it<AliTRDseedV1::kNtb; it++){
    if(!(h3=(TH3S*)gROOT->FindObject(Form("hr%s%03d_t%02d", (HasMCdata()?"MC":""), fDet, it)))){
      h3 = new TH3S(
        Form("hr%s%03d_t%02d", (HasMCdata()?"MC":""), fDet, it),
        Form(" Det[%d] t_{drift}(%2d)[bin];h*tg(#theta);tg(#phi);#Delta y[cm]", fDet, it),
        35, htgt-0.0035, htgt+0.0035, // h*tgt
        36, fExB-.18, fExB+.18,       // tgp
        60, -fDyRange, fDyRange);     // dy
    } h3->Reset();
    arr->AddAt(h3, it);
  }
  return fContainer;
}

//_______________________________________________________
void AliTRDclusterResolution::UserExec(Option_t *)
{
// Fill container histograms

  if(!(fInfo = dynamic_cast<TObjArray *>(GetInputData(1)))){
    AliError("Cluster array missing.");
    return;
  }
  if(!(fEvent = dynamic_cast<AliTRDeventInfo*>(GetInputData(2)))){
    AliError("Event Info missing.");
    return;
  }
  if(!IsCalibrated()){
    LoadCalibration();
    if(!IsCalibrated()){
      AliFatal("Loading the calibration settings failed. Check OCDB access.");
      return;
    }
    fEvent->Print();
  }
  if(!fContainer->GetEntries()) Histos();

  AliDebug(2, Form("Clusters[%d]", fInfo->GetEntriesFast()));

  Int_t det, t, np;
  Float_t x, y, z, q, dy, dydx, dzdx, cov[3], covcl[3];
  TH3S *h3(NULL); TH2I *h2(NULL);

  // define limits around ExB for which x contribution is negligible
  const Float_t kAroundZero = 3.5e-2; //(+- 2 deg)

  TObjArray *arr0 = (TObjArray*)fContainer->At(kYSys);
  TObjArray *arr1 = (TObjArray*)fContainer->At(kYRes);
  TObjArray *arr10 = (TObjArray*)arr1->At(0);
  TObjArray *arr2 = (TObjArray*)fContainer->At(kSigm);

  const AliTRDclusterInfo *cli = NULL;
  TIterator *iter=fInfo->MakeIterator();
  while((cli=dynamic_cast<AliTRDclusterInfo*>((*iter)()))){
    if((np = cli->GetNpads())>4) continue;
    cli->GetCluster(det, x, y, z, q, t, covcl);

    // select cluster according to detector region if specified
    if(fDet>=0 && fDet!=det) continue;
    if(fCol>=0 && fRow>=0){
      Int_t c,r;
      cli->GetCenterPad(c, r);
      if(TMath::Abs(fCol-c) > 5) continue;
      if(TMath::Abs(fRow-r) > 2) continue;
    }
    dy = cli->GetResolution();
    AliDebug(4, Form("det[%d] tb[%2d] q[%4.0f Log[%6.4f]] np[%d] dy[%7.2f][um] ypull[%5.2f]", det, t, q, TMath::Log(q), np, 1.e4*dy, dy/TMath::Sqrt(covcl[0])));
    
    cli->GetGlobalPosition(y, z, dydx, dzdx, &cov[0]);
    Float_t pw(cli->GetYDisplacement());

    // systematics as a function of pw and log(q)
    // only for dydx = exB + h*dzdx
    if(TMath::Abs(dydx-fExB-fH*dzdx) < kAroundZero){
      h3 = (TH3S*)arr0->At(0);
      h3->Fill(TMath::Log(q), pw, dy);
    }
    // resolution/pull as a function of pw and log(q)
    // only for dydx = 0, ExB=0
    if(TMath::Abs(fExB) < kAroundZero &&
       TMath::Abs(dydx) < kAroundZero &&
       t>=5 && t<=20 ){
      switch(np){
      case 3: // MPV np
        h3 = (TH3S*)arr10->At(t-5);
        h3->Fill(TMath::Log(q), pw, dy);
        h3 = (TH3S*)arr1->At(5);
        h3->Fill(TMath::Log(q), pw, t);
        break;
      case 2: // Min np
        h3 = (TH3S*)arr1->At(3);
        h3->Fill(TMath::Log(q), pw, dy);
        break;
      case 4: // Max np
        h3 = (TH3S*)arr1->At(4);
        h3->Fill(TMath::Log(q), pw, dy);
        break;
      }
      h3 = (TH3S*)arr1->At(1);
      h3->Fill(TMath::Log(q), pw, dy/TMath::Sqrt(covcl[0]));
    }

    // do not use problematic clusters in resolution analysis
    // TODO define limits as calibration aware (gain) !!
    //if(!AcceptableGain(fGain)) continue;
    if(q<20. || q>250.) continue;

    // systematic as a function of time bin
    // only for dydx = exB + h*dzdx and MPV q
    if(TMath::Abs(dydx-fExB-fH*dzdx) < kAroundZero){
      h2 = (TH2I*)arr0->At(1);
      h2->Fill(t, dy);
    }
    // systematic as function of tb and tgp
    // only for MPV q
    h3 = (TH3S*)arr0->At(2);
    h3->Fill(t, dydx, dy);

    // resolution/pull as a function of time bin
    // only for dydx = 0, ExB=0 and MPV q
    if(TMath::Abs(fExB) < kAroundZero &&
       TMath::Abs(dydx) < kAroundZero){
      h3 = (TH3S*)arr1->At(2);
      h3->Fill(t, dy, dy/TMath::Sqrt(covcl[0]));
    }

    // resolution as function of tb, tgp and h*tgt
    // only for MPV q
    ((TH3S*)arr2->At(t))->Fill(fH*dzdx, dydx, dy);
  }
}


//_______________________________________________________
Bool_t AliTRDclusterResolution::PostProcess()
{
// Steer processing of various cluster resolution dependences :
//
//   - process resolution dependency cluster charge
//   if(HasProcess(kYRes)) ProcessCharge();
//   - process resolution dependency on y displacement
//   if(HasProcess(kYSys)) ProcessCenterPad();
//   - process resolution dependency on drift legth and drift cell width
//   if(HasProcess(kSigm)) ProcessSigma();
//   - process systematic shift on drift legth and drift cell width
//   if(HasProcess(kMean)) ProcessMean();

  if(!fContainer) return kFALSE;
  if(!IsCalibrated()){
    AliError("Not calibrated instance.");
    return kFALSE;
  }
  TObjArray *arr = NULL;
  TTree *t=NULL;
  if(!fResults){
    TGraphErrors *g = NULL;
    fResults = new TObjArray(kNtasks);
    fResults->SetOwner();
    fResults->AddAt(arr = new TObjArray(3), kYRes);
    arr->SetOwner();
    arr->AddAt(g = new TGraphErrors(), 0);
    g->SetLineColor(kBlue); g->SetMarkerColor(kBlue);
    g->SetMarkerStyle(7); 
    arr->AddAt(g = new TGraphErrors(), 1);
    g->SetLineColor(kRed); g->SetMarkerColor(kRed);
    g->SetMarkerStyle(23); 
    arr->AddAt(g = new TGraphErrors(), 2);
    g->SetLineColor(kGreen); g->SetMarkerColor(kGreen);
    g->SetMarkerStyle(7); 

    // pad center dependence
    fResults->AddAt(arr = new TObjArray(AliTRDgeometry::kNlayer+1), kYSys);
    arr->SetOwner();
    arr->AddAt(
    t = new TTree("cent", "dy=f(y,x,ly)"), 0);
    t->Branch("ly", &fLy, "ly/B");
    t->Branch("t", &fT, "t/F");
    t->Branch("y", &fY, "y/F");
    t->Branch("m", &fR[0], "m[2]/F");
    t->Branch("s", &fR[2], "s[2]/F");
    t->Branch("pm", &fP[0], "pm[2]/F");
    t->Branch("ps", &fP[2], "ps[2]/F");
    for(Int_t il=1; il<=AliTRDgeometry::kNlayer; il++){
      arr->AddAt(g = new TGraphErrors(), il);
      g->SetLineColor(il); g->SetLineStyle(il);
      g->SetMarkerColor(il);g->SetMarkerStyle(4); 
    }


    fResults->AddAt(t = new TTree("sigm", "dy=f(dw,x,dydx)"), kSigm);
    t->Branch("t", &fT, "t/F");
    t->Branch("x", &fX, "x/F");
    t->Branch("z", &fZ, "z/F");
    t->Branch("sx", &fR[0], "sx[2]/F");
    t->Branch("sy", &fR[2], "sy[2]/F");


    fResults->AddAt(t = new TTree("mean", "dy=f(dw,x,dydx - h dzdx)"), kMean);
    t->Branch("t", &fT, "t/F");
    t->Branch("x", &fX, "x/F");
    t->Branch("z", &fZ, "z/F");
    t->Branch("dx", &fR[0], "dx[2]/F");
    t->Branch("dy", &fR[2], "dy[2]/F");
  } else {
    TObject *o = NULL;
    TIterator *iter=fResults->MakeIterator();
    while((o=((*iter)()))) o->Clear(); // maybe it is wrong but we should never reach this point
  }
  
  // process resolution dependency on charge
  if(HasProcess(kYRes)) ProcessCharge();
  
  // process resolution dependency on y displacement
  if(HasProcess(kYSys)) ProcessNormalTracks();

  // process resolution dependency on drift legth and drift cell width
  if(HasProcess(kSigm)) ProcessSigma();

  // process systematic shift on drift legth and drift cell width
  if(HasProcess(kMean)) ProcessMean();

  return kTRUE;
}

//_______________________________________________________
Bool_t AliTRDclusterResolution::LoadCalibration()
{
// Retrieve calibration parameters from OCDB, drift velocity and t0 for the detector region specified by
// a previous call to AliTRDclusterResolution::SetCalibrationRegion().

  AliCDBManager *cdb = AliCDBManager::Instance(); // check access OCDB
  if(cdb->GetRun() < 0){
    AliError("OCDB manager not properly initialized");
    return kFALSE;
  }
  // check magnetic field
  if(!TGeoGlobalMagField::Instance() || !TGeoGlobalMagField::Instance()->IsLocked()){
    AliError("Magnetic field not available.");
    return kFALSE;
  }

  AliTRDcalibDB *fCalibration  = AliTRDcalibDB::Instance();
  AliTRDCalROC  *fCalVdriftROC(fCalibration->GetVdriftROC(fDet>=0?fDet:0))
               ,*fCalT0ROC(fCalibration->GetT0ROC(fDet>=0?fDet:0));
  const AliTRDCalDet  *fCalVdriftDet = fCalibration->GetVdriftDet();
  const AliTRDCalDet  *fCalT0Det = fCalibration->GetT0Det();

  if(IsUsingCalibParam(kVdrift)){
    fVdrift = fCalVdriftDet->GetValue(fDet>=0?fDet:0);
    if(fCol>=0 && fRow>=0) fVdrift*= fCalVdriftROC->GetValue(fCol, fRow);
  }
  fExB    = AliTRDCommonParam::Instance()->GetOmegaTau(fVdrift);
  AliTRDCommonParam::Instance()->GetDiffCoeff(fDt, fDl, fVdrift);
  if(IsUsingCalibParam(kT0)){
    fT0     = fCalT0Det->GetValue(fDet>=0?fDet:0);
    if(fCol>=0 && fRow>=0) fT0 *= fCalT0ROC->GetValue(fCol, fRow);
  }
  if(IsUsingCalibParam(kGain)) fGain = (fCol>=0 && fRow>=0)?fCalibration-> GetGainFactor(fDet, fCol, fRow):fCalibration-> GetGainFactorAverage(fDet);

  SetBit(kCalibrated);

  AliInfo(Form("Calibration D[%3d] Col[%3d] Row[%2d] : \n   t0[%5.3f] vd[%5.3f] gain[%5.3f] ExB[%f] DiffT[%f] DiffL[%f]", fDet, fCol, fRow, fT0, fVdrift, fGain, fExB, fDt, fDl));

  return kTRUE;
}

//_______________________________________________________
Bool_t AliTRDclusterResolution::LoadGlobalChamberPosition()
{
// Retrieve global chamber position from alignment
// a previous call to AliTRDclusterResolution::SetCalibrationRegion() is mandatory.

  TGeoHMatrix *matrix(NULL);
  Double_t loc[] = {0., 0., 0.}, glb[] = {0., 0., 0.};
  AliTRDgeometry *geo(AliTRDinfoGen::Geometry());
  if(!(matrix= geo->GetClusterMatrix(fDet))) {
    AliFatal(Form("Could not get transformation matrix for detector %d.", fDet));
    return kFALSE;
  }
  matrix->LocalToMaster(loc, glb);
  fXch = glb[0]+ AliTRDgeometry::AnodePos()-.5*AliTRDgeometry::AmThick() - AliTRDgeometry::DrThick();
  AliTRDpadPlane *pp(geo->GetPadPlane(fDet));
  fH = TMath::Tan(pp->GetTiltingAngle()*TMath::DegToRad());
  fZch = 0.;
  if(fRow>=0){
    fZch = pp->GetRowPos(fRow)+0.5*pp->GetLengthIPad();
  }else{
    Int_t  nrows(pp->GetNrows());
    Float_t zmax(pp->GetRow0()),
            zmin(zmax - 2 * pp->GetLengthOPad()
                - (nrows-2) * pp->GetLengthIPad()
                - (nrows-1) * pp->GetRowSpacing());
    fZch = 0.5*(zmin+zmax);
  }
  
  AliInfo(Form("Global pos. D[%3d] Col[%3d] Row[%2d] : \n   x[%6.2f] z[%6.2f] h[%+6.2f].", fDet, fCol, fRow, fXch, fZch, fH));
  SetBit(kGlobal);
  return kTRUE;
}

//_______________________________________________________
void AliTRDclusterResolution::SetCalibrationRegion(Int_t det, Int_t col, Int_t row)
{
// Set calibration region in terms of detector and pad. 
// By default detector 0 mean values are considered.

  if(det>=0 && det<AliTRDgeometry::kNdet){ 
    fDet = det;
    if(col>=0 && row>=0){
      // check pad col/row for detector
      AliTRDgeometry *geo = AliTRDinfoGen::Geometry();
      AliTRDpadPlane *pp(geo->GetPadPlane(fDet));
      if(fCol>=pp->GetNcols() ||
        fRow>=pp->GetNrows()){
        AliWarning(Form("Pad coordinates col[%d] or row[%d] incorrect for det[%d].\nLimits are max col[%d] max row[%d]. Reset to default", fCol, fRow, fDet, pp->GetNcols(), pp->GetNrows()));
        fCol = -1; fRow=-1;
      } else {
        fCol = col;
        fRow = row;
      }
    }
  } else {
    AliFatal(Form("Detector index outside range [0 %d].", AliTRDgeometry::kNdet));
  }
  return;
}

//_______________________________________________________
void AliTRDclusterResolution::SetVisual()
{
  if(fCanvas) return;
  fCanvas = new TCanvas("clResCanvas", "Cluster Resolution Visualization", 10, 10, 600, 600);
}

//_______________________________________________________
void AliTRDclusterResolution::ProcessCharge()
{
// Resolution as a function of cluster charge.
//
// As described in the function ProcessCenterPad() the error parameterization for clusters for phi = a_L can be 
// written as:
// BEGIN_LATEX
// #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
// END_LATEX
// with the contribution in case of B=0 given by:
// BEGIN_LATEX
// #sigma_{y}|_{B=0} = #sigma_{diff}*Gauss(0, s_{ly}) + #delta_{#sigma}(q)
// END_LATEX
// which further can be simplified to:
// BEGIN_LATEX
// <#sigma_{y}|_{B=0}>(q) = <#sigma_{y}> + #delta_{#sigma}(q)
// <#sigma_{y}> = #int{f(q)#sigma_{y}dq}
// END_LATEX
// The results for s_y and f(q) are displayed below:
//Begin_Html
//<img src="TRD/clusterQerror.gif">
//End_Html
// The function has to extended to accomodate gain calibration scalling and errors.
//
// Author
// Alexandru Bercuci <A.Bercuci@gsi.de>



  TObjArray *arr(NULL);
  if(!(arr = (TObjArray*)fContainer->At(kYSys))) {
    AliError("Missing systematic container");
    return;
  }
  TH3S *h3s(NULL);
  if(!(h3s = (TH3S*)arr->At(0))){
    AliError("Missing systematic histo");
    return;
  }
  // PROCESS SYSTEMATIC
  Float_t tmin(6.5), tmax(20.5), tmed(0.5*(tmin+tmax));
  TGraphErrors *g[2]; TH1 *h(NULL);
  g[0] = new TGraphErrors();
  g[0]->SetMarkerStyle(24);g[0]->SetMarkerColor(kBlue);g[0]->SetLineColor(kBlue);
  g[1] = new TGraphErrors();
  g[1]->SetMarkerStyle(24);g[1]->SetMarkerColor(kRed);g[1]->SetLineColor(kRed);
  // define model for systematic shift vs pw
  TF1 fm("fm", "[0]+[1]*sin(x*[2])", -.45,.45);
  fm.SetParameter(0, 0.); fm.SetParameter(1, 1.e-2); fm.FixParameter(2, TMath::TwoPi());
  fm.SetParNames("#deltay", "#pm#delta", "2*#pi");
  h3s->GetXaxis()->SetRangeUser(tmin, tmax);
  if(!AliTRDresolution::Process((TH2*)h3s->Project3D("zy"), g))return;
  g[0]->Fit(&fm, "QR");
  if(fCanvas){
    g[0]->Draw("apl");
    fCanvas->Modified(); fCanvas->Update();
    h = g[0]->GetHistogram();
    h->SetTitle(fm.GetTitle());
    h->GetXaxis()->SetTitle("pw");h->GetXaxis()->CenterTitle();
    h->GetYaxis()->SetTitle("#Delta y[cm]");h->GetYaxis()->CenterTitle();
    if(IsSaveAs()) fCanvas->SaveAs(Form("D%03d_SysNormTrack_pw.gif", fDet));
    else gSystem->Sleep(100);
  }

  // define model for systematic shift vs tb
  TF1 fx("fx", "[0]+0.1*[1]*(x-[2])", tmin, tmax);
  fx.SetParNames("#deltay", "#deltay/t", "<t>");
  fx.FixParameter(2, tmed);
  h3s->GetXaxis()->UnZoom();
  if(!AliTRDresolution::Process((TH2*)h3s->Project3D("zx"), g)) return;
  g[0]->Fit(&fx, "Q", "", tmin, tmax);
  if(fCanvas){
    g[0]->Draw("apl");
    fCanvas->Modified(); fCanvas->Update();
    h = g[0]->GetHistogram();
    h->SetTitle(fx.GetTitle());
    h->GetXaxis()->SetTitle("t [tb]");h->GetXaxis()->CenterTitle();
    h->GetYaxis()->SetTitle("#Delta y[cm]");h->GetYaxis()->CenterTitle();
    if(IsSaveAs()) fCanvas->SaveAs(Form("D%03d_SysNormTrack_tb.gif", fDet));
    else gSystem->Sleep(100);
  }

  TH3S *h3(NULL);
  if(!(h3 = (TH3S*)fContainer->At(kYRes))) {
    AliWarning("Missing dy=f(Q) histo");
    return;
  }
  TF1 f("f", "gaus", -.5, .5);
  TAxis *ax(NULL);
  TH1 *h1(NULL);

  // compute mean error on x
  Double_t s2x = 0.; 
  for(Int_t ix=5; ix<AliTRDseedV1::kNtb; ix++){
    // retrieve error on the drift length
    s2x += AliTRDcluster::GetSX(ix);
  }
  s2x /= (AliTRDseedV1::kNtb-5); s2x *= s2x;
  //Double_t exb2 = fExB*fExB;

  arr = (TObjArray*)fResults->At(kYRes);
  TGraphErrors *gqm = (TGraphErrors*)arr->At(0);
  TGraphErrors *gqs = (TGraphErrors*)arr->At(1);
  TGraphErrors *gqp = (TGraphErrors*)arr->At(2);
  Double_t q, n = 0., entries;
  ax = h3->GetXaxis();
  for(Int_t ix=1; ix<=ax->GetNbins(); ix++){
    q = TMath::Exp(ax->GetBinCenter(ix));
    ax->SetRange(ix, ix);
    h1 = h3->Project3D("y");
    entries = h1->GetEntries();
    if(entries < 150) continue;
    h1->Fit(&f, "Q");

    // Fill sy^2 = f(q)
    Int_t ip = gqm->GetN();
    gqm->SetPoint(ip, q, 1.e4*f.GetParameter(1));
    gqm->SetPointError(ip, 0., 1.e4*f.GetParError(1));

    // correct sigma for ExB effect
    gqs->SetPoint(ip, q, 1.e4*f.GetParameter(2)/**f.GetParameter(2)-exb2*s2x)*/);
    gqs->SetPointError(ip, 0., 1.e4*f.GetParError(2)/**f.GetParameter(2)*/);

    // save probability
    n += entries;
    gqp->SetPoint(ip, q, entries);
    gqp->SetPointError(ip, 0., 0./*TMath::Sqrt(entries)*/);
  } 

  // normalize probability and get mean sy
  Double_t sm = 0., sy;
  for(Int_t ip=gqp->GetN(); ip--;){
    gqp->GetPoint(ip, q, entries);
    entries/=n;
    gqp->SetPoint(ip, q, 1.e4*entries);
    gqs->GetPoint(ip, q, sy);
    sm += entries*sy;
  }

  // error parametrization s(q) = <sy> + b(1/q-1/q0)
  TF1 fq("fq", "[0] + [1]/x", 20., 250.);
  gqs->Fit(&fq/*, "W"*/);
  printf("sm=%f [0]=%f [1]=%f\n", 1.e-4*sm, fq.GetParameter(0), fq.GetParameter(1));
  printf("  const Float_t sq0inv = %f; // [1/q0]\n", (sm-fq.GetParameter(0))/fq.GetParameter(1));
  printf("  const Float_t sqb    = %f; // [cm]\n", 1.e-4*fq.GetParameter(1));
}

//_______________________________________________________
Bool_t AliTRDclusterResolution::ProcessNormalTracks()
{
// Resolution as a function of y displacement from pad center and drift length.
//
// Since the error parameterization of cluster r-phi position can be written as (see AliTRDcluster::SetSigmaY2()):
// BEGIN_LATEX
// #sigma_{y}^{2} = (#sigma_{diff}*Gauss(0, s_{ly}) + #delta_{#sigma}(q))^{2} + tg^{2}(#alpha_{L})*#sigma_{x}^{2} + tg^{2}(#phi-#alpha_{L})*#sigma_{x}^{2}+[tg(#phi-#alpha_{L})*tg(#alpha_{L})*x]^{2}/12
// END_LATEX
// one can see that for phi = a_L one gets the following expression:
// BEGIN_LATEX
// #sigma_{y}^{2} = #sigma_{y}^{2}|_{B=0} + tg^{2}(#alpha_{L})*#sigma_{x}^{2}
// END_LATEX
// where we have explicitely marked the remaining term in case of absence of magnetic field. Thus one can use the 
// previous equation to estimate s_y for B=0 and than by comparing in magnetic field conditions one can get the s_x.
// This is a simplified method to determine the error parameterization for s_x and s_y as compared to the one 
// implemented in ProcessSigma(). For more details on cluster error parameterization please see also 
// AliTRDcluster::SetSigmaY2()
// 
// The representation of dy=f(y_cen, x_drift| layer) can be also used to estimate the systematic shift in the r-phi 
// coordinate resulting from imperfection in the cluster shape parameterization. From the expresion of the shift derived 
// in ProcessMean() with phi=exb one gets: 
// BEGIN_LATEX
// <#Delta y>= <#delta x> * (tg(#alpha_{L})-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
// <#Delta y>(y_{cen})= -h*<#delta x>(x_{drift}, q_{cl}) * dz/dx + #delta y(y_{cen}, ...)
// END_LATEX
// where all dependences are made explicit. This last expression can be used in two ways:
//   - by average on the dz/dx we can determine directly dy (the method implemented here) 
//   - by plotting as a function of dzdx one can determine both dx and dy components in an independent method.
//Begin_Html
//<img src="TRD/clusterYcorr.gif">
//End_Html
// Author
// Alexandru Bercuci <A.Bercuci@gsi.de>

  TObjArray *arr(NULL);
  TH3S *h3r(NULL), *h3t(NULL);
  if(!(arr= (TObjArray*)fContainer->At(kYRes))) {
    AliError("Missing resolution container");
    return kFALSE;
  }
  if(!(h3r = (TH3S*)arr->At(0))){
    AliError("Missing resolution pw/q histo");
    return kFALSE;
  } else if(!(Int_t)h3r->GetEntries()){
    AliError("Empty resolution pw/q histo");
    return kFALSE;
  }
  if(!(h3t = (TH3S*)arr->At(2))){
    AliError("Missing resolution t histo");
    return kFALSE;
  } else if(!(Int_t)h3t->GetEntries()){
    AliError("Empty resolution t histo");
    return kFALSE;
  }

  // local variables
  Double_t x(0.), y(0.), ex(0.), ey(0.);
  Float_t tmin(6.5), tmax(20.5), tmed(0.5*(tmin+tmax));
  TGraphErrors *g[2]; TH1 *h(NULL);
  g[0] = new TGraphErrors();
  g[0]->SetMarkerStyle(24);g[0]->SetMarkerColor(kBlue);g[0]->SetLineColor(kBlue);
  g[1] = new TGraphErrors();
  g[1]->SetMarkerStyle(24);g[1]->SetMarkerColor(kRed);g[1]->SetLineColor(kRed);

  // PROCESS RESOLUTION VS TB
  TF1 fsx("fsx", "[0]*[0]+[1]*[1]*[2]*0.1*(x-[3])", tmin, tmax);
  fsx.SetParNames("#sqrt{<#sigma^{2}(prf, q)>}(t_{med})", "D_{T}", "v_{drift}", "t_{med}");
  fsx.FixParameter(1, fDt);
  fsx.SetParameter(2, fVdrift);
  fsx.FixParameter(3, tmed);
  if(!AliTRDresolution::Process((TH2*)h3r->Project3D("yx"), g)) return kFALSE;
  for(Int_t ip(0); ip<g[1]->GetN(); ip++){
    g[1]->GetPoint(ip, x, y);ex = g[1]->GetErrorX(ip); ey = g[1]->GetErrorY(ip);
    g[1]->SetPoint(ip, x, y*y);g[1]->SetPointError(ip, ex, 2*y*ey);
  }
  g[1]->Fit(&fsx, "Q", "", tmin, tmax);
  if(fCanvas){
    g[1]->Draw("apl");
    fCanvas->Modified(); fCanvas->Update();
    h = g[1]->GetHistogram();
    h->SetTitle(fsx.GetTitle());
    h->GetXaxis()->SetTitle("t [tb]");h->GetXaxis()->CenterTitle();
    h->GetYaxis()->SetTitle("#sigma^{2} (y) [cm^{2}]");h->GetYaxis()->CenterTitle();
    if(IsSaveAs()) fCanvas->SaveAs(Form("D%03d_ResNormTrack_tb.gif", fDet));
    else gSystem->Sleep(100);
  }

  // define model for resolution vs pw
  TF1 fg("fg", "gaus", -.5, .5); fg.FixParameter(1, 0.);
  TF1 fs("fs", "[0]*[0]*exp(-1*(x/[1])**2)+[2]", -.5, .5);
  fs.SetParNames("<#sigma^{max}(q,prf)>_{q}", "#sigma(pw)", "D_{T}^{2}*<x>");
  h3r->GetXaxis()->SetRangeUser(tmin, tmax);
  if(!AliTRDresolution::Process((TH2*)h3r->Project3D("zy"), g, 200)) return kFALSE;
  for(Int_t ip(0); ip<g[1]->GetN(); ip++){
    g[1]->GetPoint(ip, x, y); ex = g[1]->GetErrorX(ip); ey = g[1]->GetErrorY(ip);
    g[1]->SetPoint(ip, x, y*y);g[1]->SetPointError(ip, ex, 2.*y*ey);
  }
  g[1]->Fit(&fg, "QR");
  fs.SetParameter(0, TMath::Sqrt(fg.GetParameter(0)));
  fs.SetParameter(1, fg.GetParameter(2));
  Float_t sdiff(fDt*fDt*fsx.GetParameter(2)*tmed*0.1);
  fs.SetParameter(2, sdiff);
  fs.SetParLimits(2, 0.1*sdiff, 1.9*sdiff);
  g[1]->Fit(&fs, "QR");
  if(fCanvas){
    g[1]->Draw("apl");
    fCanvas->Modified(); fCanvas->Update();
    h = g[1]->GetHistogram();
    h->SetTitle(fs.GetTitle());
    h->GetXaxis()->SetTitle("pw");h->GetXaxis()->CenterTitle();
    h->GetYaxis()->SetTitle("#sigma^{2} (y) [cm^{2}]");h->GetYaxis()->CenterTitle();
    if(IsSaveAs()) fCanvas->SaveAs(Form("D%03d_ResNormTrack_pw.gif", fDet));
    else gSystem->Sleep(100);
  }

  AliDebug(2, Form("<s(q,prf)>[mum] = %7.3f", 1.e4*TMath::Sqrt(fsx.Eval(0.))));
  AliDebug(2, Form("<s(q)>[mum] = %7.3f", 1.e4*TMath::Sqrt(fs.Eval(-0.5)-fs.GetParameter(2))));
  AliDebug(2, Form("<s(x)>[mum] = %7.3f(prf) %7.3f(diff)", 1.e4*TMath::Sqrt(fs.GetParameter(2)), 1.e4*TMath::Sqrt(sdiff)));

  // define model for resolution vs q
  TF1 fq("fq", "[0]*[0]*exp(-1*[1]*(x-[2])**2)+[2]", 2.5, 5.5);
  fq.SetParNames("<#sigma^{max}(q,prf)>_{prf}", "slope","mean", "D_{T}^{2}*<x>");
  if(!AliTRDresolution::Process((TH2*)h3t->Project3D("yx"), g)) return kFALSE;
  for(Int_t ip(0); ip<g[1]->GetN(); ip++){
    g[1]->GetPoint(ip, x, y); ex = g[1]->GetErrorX(ip); ey = g[1]->GetErrorY(ip);
    g[1]->SetPoint(ip, x, y*y);g[1]->SetPointError(ip, ex, 2.*y*ey);
  }
  fq.SetParameter(0, 8.e-2); fq.SetParLimits(0, 0., 1.);
  fq.SetParameter(1, 1.); //fq.SetParLimits(1, -1., 0.);
  fq.SetParameter(3, sdiff); fq.SetParLimits(3, 0.1*sdiff, 1.9*sdiff);
  g[1]->Fit(&fq, "QR");
//   AliDebug(2, Form("<sq>[mum] = %7.3f", 1.e4*TMath::Sqrt(fs.Eval(-0.5)-fs.GetParameter(2)));
//   AliDebug(2, Form("<sx>[mum] = %7.3f(prf) %7.3f(diff)", 1.e4*TMath::Sqrt(fs.Eval(-0.5)-fs.GetParameter(2)), 1.e4*TMath::Sqrt(sdiff)));
  if(fCanvas){
    g[1]->Draw("apl");
    fCanvas->Modified(); fCanvas->Update();
    h = g[1]->GetHistogram();
    h->SetTitle(fs.GetTitle());
    h->GetXaxis()->SetTitle("log(q) [a.u.]");h->GetXaxis()->CenterTitle();
    h->GetYaxis()->SetTitle("#sigma^{2} (y) [cm^{2}]");h->GetYaxis()->CenterTitle();
    if(IsSaveAs()) fCanvas->SaveAs(Form("D%03d_ResNormTrack_q.gif", fDet));
    else gSystem->Sleep(100);
  }
  return kTRUE;
}

//_______________________________________________________
void AliTRDclusterResolution::ProcessSigma()
{
// As the r-phi coordinate is the only one which is measured by the TRD detector we have to rely on it to
// estimate both the radial (x) and r-phi (y) errors. This method is based on the following assumptions. 
// The measured error in the y direction is the sum of the intrinsic contribution of the r-phi measurement
// with the contribution of the radial measurement - because x is not a parameter of Alice track model (Kalman).
// BEGIN_LATEX
// #sigma^{2}|_{y} = #sigma^{2}_{y*} + #sigma^{2}_{x*}   
// END_LATEX
// In the general case 
// BEGIN_LATEX
// #sigma^{2}_{y*} = #sigma^{2}_{y} + tg^{2}(#alpha_{L})#sigma^{2}_{x_{drift}}   
// #sigma^{2}_{x*} = tg^{2}(#phi - #alpha_{L})*(#sigma^{2}_{x_{drift}} + #sigma^{2}_{x_{0}} + tg^{2}(#alpha_{L})*x^{2}/12)
// END_LATEX
// where we have explicitely show the lorentz angle correction on y and the projection of radial component on the y
// direction through the track angle in the bending plane (phi). Also we have shown that the radial component in the
// last equation has twp terms, the drift and the misalignment (x_0). For ideal geometry or known misalignment one 
// can solve the equation
// BEGIN_LATEX
// #sigma^{2}|_{y} = tg^{2}(#phi - #alpha_{L})*(#sigma^{2}_{x} + tg^{2}(#alpha_{L})*x^{2}/12)+ [#sigma^{2}_{y} + tg^{2}(#alpha_{L})#sigma^{2}_{x}]
// END_LATEX
// by fitting a straight line:
// BEGIN_LATEX
// #sigma^{2}|_{y} = a(x_{cl}, z_{cl}) * tg^{2}(#phi - #alpha_{L}) + b(x_{cl}, z_{cl})
// END_LATEX
// the error parameterization will be given by:
// BEGIN_LATEX
// #sigma_{x} (x_{cl}, z_{cl}) = #sqrt{a(x_{cl}, z_{cl}) - tg^{2}(#alpha_{L})*x^{2}/12}
// #sigma_{y} (x_{cl}, z_{cl}) = #sqrt{b(x_{cl}, z_{cl}) - #sigma^{2}_{x} (x_{cl}, z_{cl}) * tg^{2}(#alpha_{L})}
// END_LATEX
// Below there is an example of such dependency. 
//Begin_Html
//<img src="TRD/clusterSigmaMethod.gif">
//End_Html
//
// The error parameterization obtained by this method are implemented in the functions AliTRDcluster::GetSX() and
// AliTRDcluster::GetSYdrift(). For an independent method to determine s_y as a function of drift length check the 
// function ProcessCenterPad(). One has to keep in mind that while this method return the mean s_y over the distance
// to pad center distribution the other method returns the *STANDARD* value at center=0 (maximum). To recover the 
// standard value one has to solve the obvious equation:
// BEGIN_LATEX
// #sigma_{y}^{STANDARD} = #frac{<#sigma_{y}>}{#int{s exp(s^{2}/#sigma) ds}}
// END_LATEX
// with "<s_y>" being the value calculated here and "sigma" the width of the s_y distribution calculated in 
// ProcessCenterPad().
//  
// Author
// Alexandru Bercuci <A.Bercuci@gsi.de>

  TObjArray *arr = (TObjArray*)fContainer->At(kSigm);
  if(!arr){
    AliWarning("Missing dy=f(x_d, d_w) container");
    return;
  }

  // init visualization
  TGraphErrors *ggs = NULL;
  TGraph *line = NULL;
  if(fCanvas){
    ggs = new TGraphErrors();
    line = new TGraph();
    line->SetLineColor(kRed);line->SetLineWidth(2);
  }

  // init logistic support
  TF1 f("f", "gaus", -.5, .5);
  TLinearFitter gs(1,"pol1");
  TH1 *hFrame=NULL;
  TH1D *h1 = NULL; TH3S *h3=NULL;
  TAxis *ax = NULL;
  Double_t exb2 = fExB*fExB;
  AliTRDcluster c;
  TTree *t = (TTree*)fResults->At(kSigm);
  for(Int_t ix=0; ix<AliTRDseedV1::kNtb; ix++){
    if(!(h3=(TH3S*)arr->At(ix))) continue;
    c.SetPadTime(ix);
    fX = c.GetXloc(fT0, fVdrift);
    fT = c.GetLocalTimeBin(); // ideal
    printf(" pad time[%d] local[%f]\n", ix, fT);
    for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
      ax = h3->GetXaxis();
      ax->SetRange(iz, iz);
      fZ = ax->GetBinCenter(iz);

      // reset visualization
      if(fCanvas){ 
        new(ggs) TGraphErrors();
        ggs->SetMarkerStyle(7);
      }
      gs.ClearPoints();

      for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
        ax = h3->GetYaxis();
        ax->SetRange(ip, ip); 
        Double_t tgl = ax->GetBinCenter(ip);
        // finish navigation in the HnSparse

        //if(TMath::Abs(dydx)>0.18) continue;
        Double_t tgg = (tgl-fExB)/(1.+tgl*fExB);
        Double_t tgg2 = tgg*tgg;

        h1 = (TH1D*)h3->Project3D("z");
        Int_t entries = (Int_t)h1->Integral();
        if(entries < 50) continue;
        //Adjust(&f, h1);
        h1->Fit(&f, "QN");

        Double_t s2  = f.GetParameter(2)*f.GetParameter(2);
        Double_t s2e = 2.*f.GetParameter(2)*f.GetParError(2);
        // Fill sy^2 = f(tg^2(phi-a_L))
        gs.AddPoint(&tgg2, s2, s2e);

        if(!ggs) continue;
        Int_t jp = ggs->GetN();
        ggs->SetPoint(jp, tgg2, s2);
        ggs->SetPointError(jp, 0., s2e);
      }
      // TODO here a more robust fit method has to be provided
      // for which lower boundaries on the parameters have to 
      // be imposed. Unfortunately the Minuit fit does not work 
      // for the TGraph in the case of B not 0.
      if(gs.Eval()) continue;

      fR[0] = gs.GetParameter(1) - fX*fX*exb2/12.;
      AliDebug(3, Form("  s2x+x2=%f ang=%f s2x=%f", gs.GetParameter(1), fX*fX*exb2/12., fR[0]));
      fR[0] = TMath::Max(fR[0], Float_t(4.e-4)); 

      // s^2_y  = s0^2_y + tg^2(a_L) * s^2_x
      // s0^2_y = f(D_L)*x + s_PRF^2 
      fR[2]= gs.GetParameter(0)-exb2*fR[0];
      AliDebug(3, Form("  s2y+s2x=%f s2y=%f", fR[0], fR[2]));
      fR[2] = TMath::Max(fR[2], Float_t(2.5e-5)); 
      fR[0] = TMath::Sqrt(fR[0]);
      fR[1] = .5*gs.GetParError(1)/fR[0];
      fR[2] = TMath::Sqrt(fR[2]);
      fR[3] = gs.GetParError(0)+exb2*exb2*gs.GetParError(1);
      t->Fill();
      AliDebug(2, Form("xd=%4.2f[cm] sx=%6.1f[um] sy=%5.1f[um]", fX, 1.e4*fR[0], 1.e4*fR[2]));

      if(!fCanvas) continue;
      fCanvas->cd(); fCanvas->SetLogx(); //fCanvas->SetLogy();
      if(!hFrame){ 
        fCanvas->SetMargin(0.15, 0.01, 0.1, 0.01);
        hFrame=new TH1I("hFrame", "", 100, 0., .3);
        hFrame->SetMinimum(0.);hFrame->SetMaximum(.005);
        hFrame->SetXTitle("tg^{2}(#phi-#alpha_{L})");
        hFrame->SetYTitle("#sigma^{2}y[cm^{2}]");
        hFrame->GetYaxis()->SetTitleOffset(2.);
        hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
        hFrame->Draw();
      } else hFrame->Reset();
      Double_t xx = 0., dxx=.2/50;
      for(Int_t ip=0;ip<50;ip++){ 
        line->SetPoint(ip, xx,  gs.GetParameter(0)+xx*gs.GetParameter(1)); 
        xx+=dxx;
      }
      ggs->Draw("pl"); line->Draw("l");
      fCanvas->Modified(); fCanvas->Update();
      if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessSigma_z[%5.3f]_x[%5.3f].gif", fZ, fX));
      else gSystem->Sleep(100);
    }
  }
  return;
}

//_______________________________________________________
void AliTRDclusterResolution::ProcessMean()
{
// By this method the cluster shift in r-phi and radial directions can be estimated by comparing with the MC.
// The resolution of the cluster corrected for pad tilt with respect to MC in the r-phi (measuring) plane can be 
// expressed by:
// BEGIN_LATEX
// #Delta y=w - y_{MC}(x_{cl})
// w = y_{cl}^{'} + h*(z_{MC}(x_{cl})-z_{cl})
// y_{MC}(x_{cl}) = y_{0} - dy/dx*x_{cl}
// z_{MC}(x_{cl}) = z_{0} - dz/dx*x_{cl}
// y_{cl}^{'} = y_{cl}-x_{cl}*tg(#alpha_{L})
// END_LATEX
// where x_cl is the drift length attached to a cluster, y_cl is the r-phi coordinate of the cluster measured by
// charge sharing on adjacent pads and y_0 and z_0 are MC reference points (as example the track references at 
// entrance/exit of a chamber). If we suppose that both r-phi (y) and radial (x) coordinate of the clusters are 
// affected by errors we can write
// BEGIN_LATEX
// x_{cl} = x_{cl}^{*} + #delta x 
// y_{cl} = y_{cl}^{*} + #delta y 
// END_LATEX 
// where the starred components are the corrected values. Thus by definition the following quantity
// BEGIN_LATEX
// #Delta y^{*}= w^{*} - y_{MC}(x_{cl}^{*})
// END_LATEX
// has 0 average over all dependency. Using this decomposition we can write:
// BEGIN_LATEX
// <#Delta y>=<#Delta y^{*}> + <#delta x * (dy/dx-h*dz/dx) + #delta y - #delta x * tg(#alpha_{L})>
// END_LATEX
// which can be transformed to the following linear dependence:
// BEGIN_LATEX
// <#Delta y>= <#delta x> * (dy/dx-h*dz/dx) + <#delta y - #delta x * tg(#alpha_{L})>
// END_LATEX
// if expressed as function of dy/dx-h*dz/dx. Furtheremore this expression can be plotted for various clusters
// i.e. we can explicitely introduce the diffusion (x_cl) and drift cell - anisochronity (z_cl) dependences. From 
// plotting this dependence and linear fitting it with:
// BEGIN_LATEX
// <#Delta y>= a(x_{cl}, z_{cl}) * (dy/dx-h*dz/dx) + b(x_{cl}, z_{cl})
// END_LATEX
// the systematic shifts will be given by:
// BEGIN_LATEX
// #delta x (x_{cl}, z_{cl}) = a(x_{cl}, z_{cl})
// #delta y (x_{cl}, z_{cl}) = b(x_{cl}, z_{cl}) + a(x_{cl}, z_{cl}) * tg(#alpha_{L})
// END_LATEX
// Below there is an example of such dependency. 
//Begin_Html
//<img src="TRD/clusterShiftMethod.gif">
//End_Html
//
// The occurance of the radial shift is due to the following conditions 
//   - the approximation of a constant drift velocity over the drift length (larger drift velocities close to 
//     cathode wire plane)
//   - the superposition of charge tails in the amplification region (first clusters appear to be located at the 
//     anode wire)
//   - the superposition of charge tails in the drift region (shift towards anode wire)
//   - diffusion effects which convolute with the TRF thus enlarging it
//   - approximate knowledge of the TRF (approximate measuring in test beam conditions) 
// 
// The occurance of the r-phi shift is due to the following conditions 
//   - approximate model for cluster shape (LUT)
//   - rounding-up problems
//
// The numerical results for ideal simulations for the radial and r-phi shifts are displayed below and used 
// for the cluster reconstruction (see the functions AliTRDcluster::GetXcorr() and AliTRDcluster::GetYcorr()). 
//Begin_Html
//<img src="TRD/clusterShiftX.gif">
//<img src="TRD/clusterShiftY.gif">
//End_Html
// More details can be found in the presentation given during the TRD
// software meeting at the end of 2008 and beginning of year 2009, published on indico.cern.ch.
// 
// Author 
// Alexandru Bercuci <A.Bercuci@gsi.de>


 
  TObjArray *arr = (TObjArray*)fContainer->At(kMean);
  if(!arr){
    AliWarning("Missing dy=f(x_d, d_w) container");
    return;
  }

  // init logistic support
  TF1 f("f", "gaus", -.5, .5);
  TF1 line("l", "[0]+[1]*x", -.15, .15);
  TGraphErrors *gm = new TGraphErrors();
  TH1 *hFrame=NULL;
  TH1D *h1 = NULL; TH3S *h3 =NULL;
  TAxis *ax = NULL;

  AliDebug(1, Form("Calibrate for Det[%3d] t0[%5.3f] vd[%5.3f]", fDet, fT0, fVdrift));

  AliTRDcluster c;
  TTree *t = (TTree*)fResults->At(kMean);
  for(Int_t ix=0; ix<AliTRDseedV1::kNtb; ix++){
    if(!(h3=(TH3S*)arr->At(ix))) continue;
    c.SetPadTime(ix);
    fX = c.GetXloc(fT0, fVdrift);
    fT = c.GetLocalTimeBin();
    for(Int_t iz=1; iz<=h3->GetXaxis()->GetNbins(); iz++){
      ax = h3->GetXaxis();
      ax->SetRange(iz, iz);
      fZ = ax->GetBinCenter(iz);

      // reset fitter
      new(gm) TGraphErrors();
      gm->SetMarkerStyle(7);

      for(Int_t ip=1; ip<=h3->GetYaxis()->GetNbins(); ip++){
        ax = h3->GetYaxis();
        ax->SetRange(ip, ip); 
        Double_t tgl = ax->GetBinCenter(ip);
        // finish navigation in the HnSparse

        h1 = (TH1D*)h3->Project3D("z");
        Int_t entries = (Int_t)h1->Integral();
        if(entries < 50) continue;
        //Adjust(&f, h1);
        h1->Fit(&f, "QN");

        // Fill <Dy> = f(dydx - h*dzdx)
        Int_t jp = gm->GetN();
        gm->SetPoint(jp, tgl, f.GetParameter(1));
        gm->SetPointError(jp, 0., f.GetParError(1));
      }
      if(gm->GetN()<10) continue;

      gm->Fit(&line, "QN");
      fR[0] = line.GetParameter(1); // dx
      fR[1] = line.GetParError(1);
      fR[2] = line.GetParameter(0) + fExB*fR[0]; // xs = dy - tg(a_L)*dx
      t->Fill();
      AliDebug(2, Form("tb[%02d] xd=%4.2f[cm] dx=%6.2f[um] dy=%6.2f[um]", ix, fX, 1.e4*fR[0], 1.e4*fR[2]));
      if(!fCanvas) continue;

      fCanvas->cd();
      if(!hFrame){ 
        fCanvas->SetMargin(0.1, 0.02, 0.1, 0.01);
        hFrame=new TH1I("hFrame", "", 100, -.3, .3);
        hFrame->SetMinimum(-.1);hFrame->SetMaximum(.1);
        hFrame->SetXTitle("tg#phi-htg#theta");
        hFrame->SetYTitle("#Delta y[cm]");
        hFrame->GetYaxis()->SetTitleOffset(1.5);
        hFrame->SetLineColor(1);hFrame->SetLineWidth(1);
        hFrame->Draw();
      } else hFrame->Reset();
      gm->Draw("pl"); line.Draw("same");
      fCanvas->Modified(); fCanvas->Update();
      if(IsSaveAs()) fCanvas->SaveAs(Form("Figures/ProcessMean_Z[%5.3f]_TB[%02d].gif", fZ, ix));
      else gSystem->Sleep(100);
    }
  }
}
 AliTRDclusterResolution.cxx:1
 AliTRDclusterResolution.cxx:2
 AliTRDclusterResolution.cxx:3
 AliTRDclusterResolution.cxx:4
 AliTRDclusterResolution.cxx:5
 AliTRDclusterResolution.cxx:6
 AliTRDclusterResolution.cxx:7
 AliTRDclusterResolution.cxx:8
 AliTRDclusterResolution.cxx:9
 AliTRDclusterResolution.cxx:10
 AliTRDclusterResolution.cxx:11
 AliTRDclusterResolution.cxx:12
 AliTRDclusterResolution.cxx:13
 AliTRDclusterResolution.cxx:14
 AliTRDclusterResolution.cxx:15
 AliTRDclusterResolution.cxx:16
 AliTRDclusterResolution.cxx:17
 AliTRDclusterResolution.cxx:18
 AliTRDclusterResolution.cxx:19
 AliTRDclusterResolution.cxx:20
 AliTRDclusterResolution.cxx:21
 AliTRDclusterResolution.cxx:22
 AliTRDclusterResolution.cxx:23
 AliTRDclusterResolution.cxx:24
 AliTRDclusterResolution.cxx:25
 AliTRDclusterResolution.cxx:26
 AliTRDclusterResolution.cxx:27
 AliTRDclusterResolution.cxx:28
 AliTRDclusterResolution.cxx:29
 AliTRDclusterResolution.cxx:30
 AliTRDclusterResolution.cxx:31
 AliTRDclusterResolution.cxx:32
 AliTRDclusterResolution.cxx:33
 AliTRDclusterResolution.cxx:34
 AliTRDclusterResolution.cxx:35
 AliTRDclusterResolution.cxx:36
 AliTRDclusterResolution.cxx:37
 AliTRDclusterResolution.cxx:38
 AliTRDclusterResolution.cxx:39
 AliTRDclusterResolution.cxx:40
 AliTRDclusterResolution.cxx:41
 AliTRDclusterResolution.cxx:42
 AliTRDclusterResolution.cxx:43
 AliTRDclusterResolution.cxx:44
 AliTRDclusterResolution.cxx:45
 AliTRDclusterResolution.cxx:46
 AliTRDclusterResolution.cxx:47
 AliTRDclusterResolution.cxx:48
 AliTRDclusterResolution.cxx:49
 AliTRDclusterResolution.cxx:50
 AliTRDclusterResolution.cxx:51
 AliTRDclusterResolution.cxx:52
 AliTRDclusterResolution.cxx:53
 AliTRDclusterResolution.cxx:54
 AliTRDclusterResolution.cxx:55
 AliTRDclusterResolution.cxx:56
 AliTRDclusterResolution.cxx:57
 AliTRDclusterResolution.cxx:58
 AliTRDclusterResolution.cxx:59
 AliTRDclusterResolution.cxx:60
 AliTRDclusterResolution.cxx:61
 AliTRDclusterResolution.cxx:62
 AliTRDclusterResolution.cxx:63
 AliTRDclusterResolution.cxx:64
 AliTRDclusterResolution.cxx:65
 AliTRDclusterResolution.cxx:66
 AliTRDclusterResolution.cxx:67
 AliTRDclusterResolution.cxx:68
 AliTRDclusterResolution.cxx:69
 AliTRDclusterResolution.cxx:70
 AliTRDclusterResolution.cxx:71
 AliTRDclusterResolution.cxx:72
 AliTRDclusterResolution.cxx:73
 AliTRDclusterResolution.cxx:74
 AliTRDclusterResolution.cxx:75
 AliTRDclusterResolution.cxx:76
 AliTRDclusterResolution.cxx:77
 AliTRDclusterResolution.cxx:78
 AliTRDclusterResolution.cxx:79
 AliTRDclusterResolution.cxx:80
 AliTRDclusterResolution.cxx:81
 AliTRDclusterResolution.cxx:82
 AliTRDclusterResolution.cxx:83
 AliTRDclusterResolution.cxx:84
 AliTRDclusterResolution.cxx:85
 AliTRDclusterResolution.cxx:86
 AliTRDclusterResolution.cxx:87
 AliTRDclusterResolution.cxx:88
 AliTRDclusterResolution.cxx:89
 AliTRDclusterResolution.cxx:90
 AliTRDclusterResolution.cxx:91
 AliTRDclusterResolution.cxx:92
 AliTRDclusterResolution.cxx:93
 AliTRDclusterResolution.cxx:94
 AliTRDclusterResolution.cxx:95
 AliTRDclusterResolution.cxx:96
 AliTRDclusterResolution.cxx:97
 AliTRDclusterResolution.cxx:98
 AliTRDclusterResolution.cxx:99
 AliTRDclusterResolution.cxx:100
 AliTRDclusterResolution.cxx:101
 AliTRDclusterResolution.cxx:102
 AliTRDclusterResolution.cxx:103
 AliTRDclusterResolution.cxx:104
 AliTRDclusterResolution.cxx:105
 AliTRDclusterResolution.cxx:106
 AliTRDclusterResolution.cxx:107
 AliTRDclusterResolution.cxx:108
 AliTRDclusterResolution.cxx:109
 AliTRDclusterResolution.cxx:110
 AliTRDclusterResolution.cxx:111
 AliTRDclusterResolution.cxx:112
 AliTRDclusterResolution.cxx:113
 AliTRDclusterResolution.cxx:114
 AliTRDclusterResolution.cxx:115
 AliTRDclusterResolution.cxx:116
 AliTRDclusterResolution.cxx:117
 AliTRDclusterResolution.cxx:118
 AliTRDclusterResolution.cxx:119
 AliTRDclusterResolution.cxx:120
 AliTRDclusterResolution.cxx:121
 AliTRDclusterResolution.cxx:122
 AliTRDclusterResolution.cxx:123
 AliTRDclusterResolution.cxx:124
 AliTRDclusterResolution.cxx:125
 AliTRDclusterResolution.cxx:126
 AliTRDclusterResolution.cxx:127
 AliTRDclusterResolution.cxx:128
 AliTRDclusterResolution.cxx:129
 AliTRDclusterResolution.cxx:130
 AliTRDclusterResolution.cxx:131
 AliTRDclusterResolution.cxx:132
 AliTRDclusterResolution.cxx:133
 AliTRDclusterResolution.cxx:134
 AliTRDclusterResolution.cxx:135
 AliTRDclusterResolution.cxx:136
 AliTRDclusterResolution.cxx:137
 AliTRDclusterResolution.cxx:138
 AliTRDclusterResolution.cxx:139
 AliTRDclusterResolution.cxx:140
 AliTRDclusterResolution.cxx:141
 AliTRDclusterResolution.cxx:142
 AliTRDclusterResolution.cxx:143
 AliTRDclusterResolution.cxx:144
 AliTRDclusterResolution.cxx:145
 AliTRDclusterResolution.cxx:146
 AliTRDclusterResolution.cxx:147
 AliTRDclusterResolution.cxx:148
 AliTRDclusterResolution.cxx:149
 AliTRDclusterResolution.cxx:150
 AliTRDclusterResolution.cxx:151
 AliTRDclusterResolution.cxx:152
 AliTRDclusterResolution.cxx:153
 AliTRDclusterResolution.cxx:154
 AliTRDclusterResolution.cxx:155
 AliTRDclusterResolution.cxx:156
 AliTRDclusterResolution.cxx:157
 AliTRDclusterResolution.cxx:158
 AliTRDclusterResolution.cxx:159
 AliTRDclusterResolution.cxx:160
 AliTRDclusterResolution.cxx:161
 AliTRDclusterResolution.cxx:162
 AliTRDclusterResolution.cxx:163
 AliTRDclusterResolution.cxx:164
 AliTRDclusterResolution.cxx:165
 AliTRDclusterResolution.cxx:166
 AliTRDclusterResolution.cxx:167
 AliTRDclusterResolution.cxx:168
 AliTRDclusterResolution.cxx:169
 AliTRDclusterResolution.cxx:170
 AliTRDclusterResolution.cxx:171
 AliTRDclusterResolution.cxx:172
 AliTRDclusterResolution.cxx:173
 AliTRDclusterResolution.cxx:174
 AliTRDclusterResolution.cxx:175
 AliTRDclusterResolution.cxx:176
 AliTRDclusterResolution.cxx:177
 AliTRDclusterResolution.cxx:178
 AliTRDclusterResolution.cxx:179
 AliTRDclusterResolution.cxx:180
 AliTRDclusterResolution.cxx:181
 AliTRDclusterResolution.cxx:182
 AliTRDclusterResolution.cxx:183
 AliTRDclusterResolution.cxx:184
 AliTRDclusterResolution.cxx:185
 AliTRDclusterResolution.cxx:186
 AliTRDclusterResolution.cxx:187
 AliTRDclusterResolution.cxx:188
 AliTRDclusterResolution.cxx:189
 AliTRDclusterResolution.cxx:190
 AliTRDclusterResolution.cxx:191
 AliTRDclusterResolution.cxx:192
 AliTRDclusterResolution.cxx:193
 AliTRDclusterResolution.cxx:194
 AliTRDclusterResolution.cxx:195
 AliTRDclusterResolution.cxx:196
 AliTRDclusterResolution.cxx:197
 AliTRDclusterResolution.cxx:198
 AliTRDclusterResolution.cxx:199
 AliTRDclusterResolution.cxx:200
 AliTRDclusterResolution.cxx:201
 AliTRDclusterResolution.cxx:202
 AliTRDclusterResolution.cxx:203
 AliTRDclusterResolution.cxx:204
 AliTRDclusterResolution.cxx:205
 AliTRDclusterResolution.cxx:206
 AliTRDclusterResolution.cxx:207
 AliTRDclusterResolution.cxx:208
 AliTRDclusterResolution.cxx:209
 AliTRDclusterResolution.cxx:210
 AliTRDclusterResolution.cxx:211
 AliTRDclusterResolution.cxx:212
 AliTRDclusterResolution.cxx:213
 AliTRDclusterResolution.cxx:214
 AliTRDclusterResolution.cxx:215
 AliTRDclusterResolution.cxx:216
 AliTRDclusterResolution.cxx:217
 AliTRDclusterResolution.cxx:218
 AliTRDclusterResolution.cxx:219
 AliTRDclusterResolution.cxx:220
 AliTRDclusterResolution.cxx:221
 AliTRDclusterResolution.cxx:222
 AliTRDclusterResolution.cxx:223
 AliTRDclusterResolution.cxx:224
 AliTRDclusterResolution.cxx:225
 AliTRDclusterResolution.cxx:226
 AliTRDclusterResolution.cxx:227
 AliTRDclusterResolution.cxx:228
 AliTRDclusterResolution.cxx:229
 AliTRDclusterResolution.cxx:230
 AliTRDclusterResolution.cxx:231
 AliTRDclusterResolution.cxx:232
 AliTRDclusterResolution.cxx:233
 AliTRDclusterResolution.cxx:234
 AliTRDclusterResolution.cxx:235
 AliTRDclusterResolution.cxx:236
 AliTRDclusterResolution.cxx:237
 AliTRDclusterResolution.cxx:238
 AliTRDclusterResolution.cxx:239
 AliTRDclusterResolution.cxx:240
 AliTRDclusterResolution.cxx:241
 AliTRDclusterResolution.cxx:242
 AliTRDclusterResolution.cxx:243
 AliTRDclusterResolution.cxx:244
 AliTRDclusterResolution.cxx:245
 AliTRDclusterResolution.cxx:246
 AliTRDclusterResolution.cxx:247
 AliTRDclusterResolution.cxx:248
 AliTRDclusterResolution.cxx:249
 AliTRDclusterResolution.cxx:250
 AliTRDclusterResolution.cxx:251
 AliTRDclusterResolution.cxx:252
 AliTRDclusterResolution.cxx:253
 AliTRDclusterResolution.cxx:254
 AliTRDclusterResolution.cxx:255
 AliTRDclusterResolution.cxx:256
 AliTRDclusterResolution.cxx:257
 AliTRDclusterResolution.cxx:258
 AliTRDclusterResolution.cxx:259
 AliTRDclusterResolution.cxx:260
 AliTRDclusterResolution.cxx:261
 AliTRDclusterResolution.cxx:262
 AliTRDclusterResolution.cxx:263
 AliTRDclusterResolution.cxx:264
 AliTRDclusterResolution.cxx:265
 AliTRDclusterResolution.cxx:266
 AliTRDclusterResolution.cxx:267
 AliTRDclusterResolution.cxx:268
 AliTRDclusterResolution.cxx:269
 AliTRDclusterResolution.cxx:270
 AliTRDclusterResolution.cxx:271
 AliTRDclusterResolution.cxx:272
 AliTRDclusterResolution.cxx:273
 AliTRDclusterResolution.cxx:274
 AliTRDclusterResolution.cxx:275
 AliTRDclusterResolution.cxx:276
 AliTRDclusterResolution.cxx:277
 AliTRDclusterResolution.cxx:278
 AliTRDclusterResolution.cxx:279
 AliTRDclusterResolution.cxx:280
 AliTRDclusterResolution.cxx:281
 AliTRDclusterResolution.cxx:282
 AliTRDclusterResolution.cxx:283
 AliTRDclusterResolution.cxx:284
 AliTRDclusterResolution.cxx:285
 AliTRDclusterResolution.cxx:286
 AliTRDclusterResolution.cxx:287
 AliTRDclusterResolution.cxx:288
 AliTRDclusterResolution.cxx:289
 AliTRDclusterResolution.cxx:290
 AliTRDclusterResolution.cxx:291
 AliTRDclusterResolution.cxx:292
 AliTRDclusterResolution.cxx:293
 AliTRDclusterResolution.cxx:294
 AliTRDclusterResolution.cxx:295
 AliTRDclusterResolution.cxx:296
 AliTRDclusterResolution.cxx:297
 AliTRDclusterResolution.cxx:298
 AliTRDclusterResolution.cxx:299
 AliTRDclusterResolution.cxx:300
 AliTRDclusterResolution.cxx:301
 AliTRDclusterResolution.cxx:302
 AliTRDclusterResolution.cxx:303
 AliTRDclusterResolution.cxx:304
 AliTRDclusterResolution.cxx:305
 AliTRDclusterResolution.cxx:306
 AliTRDclusterResolution.cxx:307
 AliTRDclusterResolution.cxx:308
 AliTRDclusterResolution.cxx:309
 AliTRDclusterResolution.cxx:310
 AliTRDclusterResolution.cxx:311
 AliTRDclusterResolution.cxx:312
 AliTRDclusterResolution.cxx:313
 AliTRDclusterResolution.cxx:314
 AliTRDclusterResolution.cxx:315
 AliTRDclusterResolution.cxx:316
 AliTRDclusterResolution.cxx:317
 AliTRDclusterResolution.cxx:318
 AliTRDclusterResolution.cxx:319
 AliTRDclusterResolution.cxx:320
 AliTRDclusterResolution.cxx:321
 AliTRDclusterResolution.cxx:322
 AliTRDclusterResolution.cxx:323
 AliTRDclusterResolution.cxx:324
 AliTRDclusterResolution.cxx:325
 AliTRDclusterResolution.cxx:326
 AliTRDclusterResolution.cxx:327
 AliTRDclusterResolution.cxx:328
 AliTRDclusterResolution.cxx:329
 AliTRDclusterResolution.cxx:330
 AliTRDclusterResolution.cxx:331
 AliTRDclusterResolution.cxx:332
 AliTRDclusterResolution.cxx:333
 AliTRDclusterResolution.cxx:334
 AliTRDclusterResolution.cxx:335
 AliTRDclusterResolution.cxx:336
 AliTRDclusterResolution.cxx:337
 AliTRDclusterResolution.cxx:338
 AliTRDclusterResolution.cxx:339
 AliTRDclusterResolution.cxx:340
 AliTRDclusterResolution.cxx:341
 AliTRDclusterResolution.cxx:342
 AliTRDclusterResolution.cxx:343
 AliTRDclusterResolution.cxx:344
 AliTRDclusterResolution.cxx:345
 AliTRDclusterResolution.cxx:346
 AliTRDclusterResolution.cxx:347
 AliTRDclusterResolution.cxx:348
 AliTRDclusterResolution.cxx:349
 AliTRDclusterResolution.cxx:350
 AliTRDclusterResolution.cxx:351
 AliTRDclusterResolution.cxx:352
 AliTRDclusterResolution.cxx:353
 AliTRDclusterResolution.cxx:354
 AliTRDclusterResolution.cxx:355
 AliTRDclusterResolution.cxx:356
 AliTRDclusterResolution.cxx:357
 AliTRDclusterResolution.cxx:358
 AliTRDclusterResolution.cxx:359
 AliTRDclusterResolution.cxx:360
 AliTRDclusterResolution.cxx:361
 AliTRDclusterResolution.cxx:362
 AliTRDclusterResolution.cxx:363
 AliTRDclusterResolution.cxx:364
 AliTRDclusterResolution.cxx:365
 AliTRDclusterResolution.cxx:366
 AliTRDclusterResolution.cxx:367
 AliTRDclusterResolution.cxx:368
 AliTRDclusterResolution.cxx:369
 AliTRDclusterResolution.cxx:370
 AliTRDclusterResolution.cxx:371
 AliTRDclusterResolution.cxx:372
 AliTRDclusterResolution.cxx:373
 AliTRDclusterResolution.cxx:374
 AliTRDclusterResolution.cxx:375
 AliTRDclusterResolution.cxx:376
 AliTRDclusterResolution.cxx:377
 AliTRDclusterResolution.cxx:378
 AliTRDclusterResolution.cxx:379
 AliTRDclusterResolution.cxx:380
 AliTRDclusterResolution.cxx:381
 AliTRDclusterResolution.cxx:382
 AliTRDclusterResolution.cxx:383
 AliTRDclusterResolution.cxx:384
 AliTRDclusterResolution.cxx:385
 AliTRDclusterResolution.cxx:386
 AliTRDclusterResolution.cxx:387
 AliTRDclusterResolution.cxx:388
 AliTRDclusterResolution.cxx:389
 AliTRDclusterResolution.cxx:390
 AliTRDclusterResolution.cxx:391
 AliTRDclusterResolution.cxx:392
 AliTRDclusterResolution.cxx:393
 AliTRDclusterResolution.cxx:394
 AliTRDclusterResolution.cxx:395
 AliTRDclusterResolution.cxx:396
 AliTRDclusterResolution.cxx:397
 AliTRDclusterResolution.cxx:398
 AliTRDclusterResolution.cxx:399
 AliTRDclusterResolution.cxx:400
 AliTRDclusterResolution.cxx:401
 AliTRDclusterResolution.cxx:402
 AliTRDclusterResolution.cxx:403
 AliTRDclusterResolution.cxx:404
 AliTRDclusterResolution.cxx:405
 AliTRDclusterResolution.cxx:406
 AliTRDclusterResolution.cxx:407
 AliTRDclusterResolution.cxx:408
 AliTRDclusterResolution.cxx:409
 AliTRDclusterResolution.cxx:410
 AliTRDclusterResolution.cxx:411
 AliTRDclusterResolution.cxx:412
 AliTRDclusterResolution.cxx:413
 AliTRDclusterResolution.cxx:414
 AliTRDclusterResolution.cxx:415
 AliTRDclusterResolution.cxx:416
 AliTRDclusterResolution.cxx:417
 AliTRDclusterResolution.cxx:418
 AliTRDclusterResolution.cxx:419
 AliTRDclusterResolution.cxx:420
 AliTRDclusterResolution.cxx:421
 AliTRDclusterResolution.cxx:422
 AliTRDclusterResolution.cxx:423
 AliTRDclusterResolution.cxx:424
 AliTRDclusterResolution.cxx:425
 AliTRDclusterResolution.cxx:426
 AliTRDclusterResolution.cxx:427
 AliTRDclusterResolution.cxx:428
 AliTRDclusterResolution.cxx:429
 AliTRDclusterResolution.cxx:430
 AliTRDclusterResolution.cxx:431
 AliTRDclusterResolution.cxx:432
 AliTRDclusterResolution.cxx:433
 AliTRDclusterResolution.cxx:434
 AliTRDclusterResolution.cxx:435
 AliTRDclusterResolution.cxx:436
 AliTRDclusterResolution.cxx:437
 AliTRDclusterResolution.cxx:438
 AliTRDclusterResolution.cxx:439
 AliTRDclusterResolution.cxx:440
 AliTRDclusterResolution.cxx:441
 AliTRDclusterResolution.cxx:442
 AliTRDclusterResolution.cxx:443
 AliTRDclusterResolution.cxx:444
 AliTRDclusterResolution.cxx:445
 AliTRDclusterResolution.cxx:446
 AliTRDclusterResolution.cxx:447
 AliTRDclusterResolution.cxx:448
 AliTRDclusterResolution.cxx:449
 AliTRDclusterResolution.cxx:450
 AliTRDclusterResolution.cxx:451
 AliTRDclusterResolution.cxx:452
 AliTRDclusterResolution.cxx:453
 AliTRDclusterResolution.cxx:454
 AliTRDclusterResolution.cxx:455
 AliTRDclusterResolution.cxx:456
 AliTRDclusterResolution.cxx:457
 AliTRDclusterResolution.cxx:458
 AliTRDclusterResolution.cxx:459
 AliTRDclusterResolution.cxx:460
 AliTRDclusterResolution.cxx:461
 AliTRDclusterResolution.cxx:462
 AliTRDclusterResolution.cxx:463
 AliTRDclusterResolution.cxx:464
 AliTRDclusterResolution.cxx:465
 AliTRDclusterResolution.cxx:466
 AliTRDclusterResolution.cxx:467
 AliTRDclusterResolution.cxx:468
 AliTRDclusterResolution.cxx:469
 AliTRDclusterResolution.cxx:470
 AliTRDclusterResolution.cxx:471
 AliTRDclusterResolution.cxx:472
 AliTRDclusterResolution.cxx:473
 AliTRDclusterResolution.cxx:474
 AliTRDclusterResolution.cxx:475
 AliTRDclusterResolution.cxx:476
 AliTRDclusterResolution.cxx:477
 AliTRDclusterResolution.cxx:478
 AliTRDclusterResolution.cxx:479
 AliTRDclusterResolution.cxx:480
 AliTRDclusterResolution.cxx:481
 AliTRDclusterResolution.cxx:482
 AliTRDclusterResolution.cxx:483
 AliTRDclusterResolution.cxx:484
 AliTRDclusterResolution.cxx:485
 AliTRDclusterResolution.cxx:486
 AliTRDclusterResolution.cxx:487
 AliTRDclusterResolution.cxx:488
 AliTRDclusterResolution.cxx:489
 AliTRDclusterResolution.cxx:490
 AliTRDclusterResolution.cxx:491
 AliTRDclusterResolution.cxx:492
 AliTRDclusterResolution.cxx:493
 AliTRDclusterResolution.cxx:494
 AliTRDclusterResolution.cxx:495
 AliTRDclusterResolution.cxx:496
 AliTRDclusterResolution.cxx:497
 AliTRDclusterResolution.cxx:498
 AliTRDclusterResolution.cxx:499
 AliTRDclusterResolution.cxx:500
 AliTRDclusterResolution.cxx:501
 AliTRDclusterResolution.cxx:502
 AliTRDclusterResolution.cxx:503
 AliTRDclusterResolution.cxx:504
 AliTRDclusterResolution.cxx:505
 AliTRDclusterResolution.cxx:506
 AliTRDclusterResolution.cxx:507
 AliTRDclusterResolution.cxx:508
 AliTRDclusterResolution.cxx:509
 AliTRDclusterResolution.cxx:510
 AliTRDclusterResolution.cxx:511
 AliTRDclusterResolution.cxx:512
 AliTRDclusterResolution.cxx:513
 AliTRDclusterResolution.cxx:514
 AliTRDclusterResolution.cxx:515
 AliTRDclusterResolution.cxx:516
 AliTRDclusterResolution.cxx:517
 AliTRDclusterResolution.cxx:518
 AliTRDclusterResolution.cxx:519
 AliTRDclusterResolution.cxx:520
 AliTRDclusterResolution.cxx:521
 AliTRDclusterResolution.cxx:522
 AliTRDclusterResolution.cxx:523
 AliTRDclusterResolution.cxx:524
 AliTRDclusterResolution.cxx:525
 AliTRDclusterResolution.cxx:526
 AliTRDclusterResolution.cxx:527
 AliTRDclusterResolution.cxx:528
 AliTRDclusterResolution.cxx:529
 AliTRDclusterResolution.cxx:530
 AliTRDclusterResolution.cxx:531
 AliTRDclusterResolution.cxx:532
 AliTRDclusterResolution.cxx:533
 AliTRDclusterResolution.cxx:534
 AliTRDclusterResolution.cxx:535
 AliTRDclusterResolution.cxx:536
 AliTRDclusterResolution.cxx:537
 AliTRDclusterResolution.cxx:538
 AliTRDclusterResolution.cxx:539
 AliTRDclusterResolution.cxx:540
 AliTRDclusterResolution.cxx:541
 AliTRDclusterResolution.cxx:542
 AliTRDclusterResolution.cxx:543
 AliTRDclusterResolution.cxx:544
 AliTRDclusterResolution.cxx:545
 AliTRDclusterResolution.cxx:546
 AliTRDclusterResolution.cxx:547
 AliTRDclusterResolution.cxx:548
 AliTRDclusterResolution.cxx:549
 AliTRDclusterResolution.cxx:550
 AliTRDclusterResolution.cxx:551
 AliTRDclusterResolution.cxx:552
 AliTRDclusterResolution.cxx:553
 AliTRDclusterResolution.cxx:554
 AliTRDclusterResolution.cxx:555
 AliTRDclusterResolution.cxx:556
 AliTRDclusterResolution.cxx:557
 AliTRDclusterResolution.cxx:558
 AliTRDclusterResolution.cxx:559
 AliTRDclusterResolution.cxx:560
 AliTRDclusterResolution.cxx:561
 AliTRDclusterResolution.cxx:562
 AliTRDclusterResolution.cxx:563
 AliTRDclusterResolution.cxx:564
 AliTRDclusterResolution.cxx:565
 AliTRDclusterResolution.cxx:566
 AliTRDclusterResolution.cxx:567
 AliTRDclusterResolution.cxx:568
 AliTRDclusterResolution.cxx:569
 AliTRDclusterResolution.cxx:570
 AliTRDclusterResolution.cxx:571
 AliTRDclusterResolution.cxx:572
 AliTRDclusterResolution.cxx:573
 AliTRDclusterResolution.cxx:574
 AliTRDclusterResolution.cxx:575
 AliTRDclusterResolution.cxx:576
 AliTRDclusterResolution.cxx:577
 AliTRDclusterResolution.cxx:578
 AliTRDclusterResolution.cxx:579
 AliTRDclusterResolution.cxx:580
 AliTRDclusterResolution.cxx:581
 AliTRDclusterResolution.cxx:582
 AliTRDclusterResolution.cxx:583
 AliTRDclusterResolution.cxx:584
 AliTRDclusterResolution.cxx:585
 AliTRDclusterResolution.cxx:586
 AliTRDclusterResolution.cxx:587
 AliTRDclusterResolution.cxx:588
 AliTRDclusterResolution.cxx:589
 AliTRDclusterResolution.cxx:590
 AliTRDclusterResolution.cxx:591
 AliTRDclusterResolution.cxx:592
 AliTRDclusterResolution.cxx:593
 AliTRDclusterResolution.cxx:594
 AliTRDclusterResolution.cxx:595
 AliTRDclusterResolution.cxx:596
 AliTRDclusterResolution.cxx:597
 AliTRDclusterResolution.cxx:598
 AliTRDclusterResolution.cxx:599
 AliTRDclusterResolution.cxx:600
 AliTRDclusterResolution.cxx:601
 AliTRDclusterResolution.cxx:602
 AliTRDclusterResolution.cxx:603
 AliTRDclusterResolution.cxx:604
 AliTRDclusterResolution.cxx:605
 AliTRDclusterResolution.cxx:606
 AliTRDclusterResolution.cxx:607
 AliTRDclusterResolution.cxx:608
 AliTRDclusterResolution.cxx:609
 AliTRDclusterResolution.cxx:610
 AliTRDclusterResolution.cxx:611
 AliTRDclusterResolution.cxx:612
 AliTRDclusterResolution.cxx:613
 AliTRDclusterResolution.cxx:614
 AliTRDclusterResolution.cxx:615
 AliTRDclusterResolution.cxx:616
 AliTRDclusterResolution.cxx:617
 AliTRDclusterResolution.cxx:618
 AliTRDclusterResolution.cxx:619
 AliTRDclusterResolution.cxx:620
 AliTRDclusterResolution.cxx:621
 AliTRDclusterResolution.cxx:622
 AliTRDclusterResolution.cxx:623
 AliTRDclusterResolution.cxx:624
 AliTRDclusterResolution.cxx:625
 AliTRDclusterResolution.cxx:626
 AliTRDclusterResolution.cxx:627
 AliTRDclusterResolution.cxx:628
 AliTRDclusterResolution.cxx:629
 AliTRDclusterResolution.cxx:630
 AliTRDclusterResolution.cxx:631
 AliTRDclusterResolution.cxx:632
 AliTRDclusterResolution.cxx:633
 AliTRDclusterResolution.cxx:634
 AliTRDclusterResolution.cxx:635
 AliTRDclusterResolution.cxx:636
 AliTRDclusterResolution.cxx:637
 AliTRDclusterResolution.cxx:638
 AliTRDclusterResolution.cxx:639
 AliTRDclusterResolution.cxx:640
 AliTRDclusterResolution.cxx:641
 AliTRDclusterResolution.cxx:642
 AliTRDclusterResolution.cxx:643
 AliTRDclusterResolution.cxx:644
 AliTRDclusterResolution.cxx:645
 AliTRDclusterResolution.cxx:646
 AliTRDclusterResolution.cxx:647
 AliTRDclusterResolution.cxx:648
 AliTRDclusterResolution.cxx:649
 AliTRDclusterResolution.cxx:650
 AliTRDclusterResolution.cxx:651
 AliTRDclusterResolution.cxx:652
 AliTRDclusterResolution.cxx:653
 AliTRDclusterResolution.cxx:654
 AliTRDclusterResolution.cxx:655
 AliTRDclusterResolution.cxx:656
 AliTRDclusterResolution.cxx:657
 AliTRDclusterResolution.cxx:658
 AliTRDclusterResolution.cxx:659
 AliTRDclusterResolution.cxx:660
 AliTRDclusterResolution.cxx:661
 AliTRDclusterResolution.cxx:662
 AliTRDclusterResolution.cxx:663
 AliTRDclusterResolution.cxx:664
 AliTRDclusterResolution.cxx:665
 AliTRDclusterResolution.cxx:666
 AliTRDclusterResolution.cxx:667
 AliTRDclusterResolution.cxx:668
 AliTRDclusterResolution.cxx:669
 AliTRDclusterResolution.cxx:670
 AliTRDclusterResolution.cxx:671
 AliTRDclusterResolution.cxx:672
 AliTRDclusterResolution.cxx:673
 AliTRDclusterResolution.cxx:674
 AliTRDclusterResolution.cxx:675
 AliTRDclusterResolution.cxx:676
 AliTRDclusterResolution.cxx:677
 AliTRDclusterResolution.cxx:678
 AliTRDclusterResolution.cxx:679
 AliTRDclusterResolution.cxx:680
 AliTRDclusterResolution.cxx:681
 AliTRDclusterResolution.cxx:682
 AliTRDclusterResolution.cxx:683
 AliTRDclusterResolution.cxx:684
 AliTRDclusterResolution.cxx:685
 AliTRDclusterResolution.cxx:686
 AliTRDclusterResolution.cxx:687
 AliTRDclusterResolution.cxx:688
 AliTRDclusterResolution.cxx:689
 AliTRDclusterResolution.cxx:690
 AliTRDclusterResolution.cxx:691
 AliTRDclusterResolution.cxx:692
 AliTRDclusterResolution.cxx:693
 AliTRDclusterResolution.cxx:694
 AliTRDclusterResolution.cxx:695
 AliTRDclusterResolution.cxx:696
 AliTRDclusterResolution.cxx:697
 AliTRDclusterResolution.cxx:698
 AliTRDclusterResolution.cxx:699
 AliTRDclusterResolution.cxx:700
 AliTRDclusterResolution.cxx:701
 AliTRDclusterResolution.cxx:702
 AliTRDclusterResolution.cxx:703
 AliTRDclusterResolution.cxx:704
 AliTRDclusterResolution.cxx:705
 AliTRDclusterResolution.cxx:706
 AliTRDclusterResolution.cxx:707
 AliTRDclusterResolution.cxx:708
 AliTRDclusterResolution.cxx:709
 AliTRDclusterResolution.cxx:710
 AliTRDclusterResolution.cxx:711
 AliTRDclusterResolution.cxx:712
 AliTRDclusterResolution.cxx:713
 AliTRDclusterResolution.cxx:714
 AliTRDclusterResolution.cxx:715
 AliTRDclusterResolution.cxx:716
 AliTRDclusterResolution.cxx:717
 AliTRDclusterResolution.cxx:718
 AliTRDclusterResolution.cxx:719
 AliTRDclusterResolution.cxx:720
 AliTRDclusterResolution.cxx:721
 AliTRDclusterResolution.cxx:722
 AliTRDclusterResolution.cxx:723
 AliTRDclusterResolution.cxx:724
 AliTRDclusterResolution.cxx:725
 AliTRDclusterResolution.cxx:726
 AliTRDclusterResolution.cxx:727
 AliTRDclusterResolution.cxx:728
 AliTRDclusterResolution.cxx:729
 AliTRDclusterResolution.cxx:730
 AliTRDclusterResolution.cxx:731
 AliTRDclusterResolution.cxx:732
 AliTRDclusterResolution.cxx:733
 AliTRDclusterResolution.cxx:734
 AliTRDclusterResolution.cxx:735
 AliTRDclusterResolution.cxx:736
 AliTRDclusterResolution.cxx:737
 AliTRDclusterResolution.cxx:738
 AliTRDclusterResolution.cxx:739
 AliTRDclusterResolution.cxx:740
 AliTRDclusterResolution.cxx:741
 AliTRDclusterResolution.cxx:742
 AliTRDclusterResolution.cxx:743
 AliTRDclusterResolution.cxx:744
 AliTRDclusterResolution.cxx:745
 AliTRDclusterResolution.cxx:746
 AliTRDclusterResolution.cxx:747
 AliTRDclusterResolution.cxx:748
 AliTRDclusterResolution.cxx:749
 AliTRDclusterResolution.cxx:750
 AliTRDclusterResolution.cxx:751
 AliTRDclusterResolution.cxx:752
 AliTRDclusterResolution.cxx:753
 AliTRDclusterResolution.cxx:754
 AliTRDclusterResolution.cxx:755
 AliTRDclusterResolution.cxx:756
 AliTRDclusterResolution.cxx:757
 AliTRDclusterResolution.cxx:758
 AliTRDclusterResolution.cxx:759
 AliTRDclusterResolution.cxx:760
 AliTRDclusterResolution.cxx:761
 AliTRDclusterResolution.cxx:762
 AliTRDclusterResolution.cxx:763
 AliTRDclusterResolution.cxx:764
 AliTRDclusterResolution.cxx:765
 AliTRDclusterResolution.cxx:766
 AliTRDclusterResolution.cxx:767
 AliTRDclusterResolution.cxx:768
 AliTRDclusterResolution.cxx:769
 AliTRDclusterResolution.cxx:770
 AliTRDclusterResolution.cxx:771
 AliTRDclusterResolution.cxx:772
 AliTRDclusterResolution.cxx:773
 AliTRDclusterResolution.cxx:774
 AliTRDclusterResolution.cxx:775
 AliTRDclusterResolution.cxx:776
 AliTRDclusterResolution.cxx:777
 AliTRDclusterResolution.cxx:778
 AliTRDclusterResolution.cxx:779
 AliTRDclusterResolution.cxx:780
 AliTRDclusterResolution.cxx:781
 AliTRDclusterResolution.cxx:782
 AliTRDclusterResolution.cxx:783
 AliTRDclusterResolution.cxx:784
 AliTRDclusterResolution.cxx:785
 AliTRDclusterResolution.cxx:786
 AliTRDclusterResolution.cxx:787
 AliTRDclusterResolution.cxx:788
 AliTRDclusterResolution.cxx:789
 AliTRDclusterResolution.cxx:790
 AliTRDclusterResolution.cxx:791
 AliTRDclusterResolution.cxx:792
 AliTRDclusterResolution.cxx:793
 AliTRDclusterResolution.cxx:794
 AliTRDclusterResolution.cxx:795
 AliTRDclusterResolution.cxx:796
 AliTRDclusterResolution.cxx:797
 AliTRDclusterResolution.cxx:798
 AliTRDclusterResolution.cxx:799
 AliTRDclusterResolution.cxx:800
 AliTRDclusterResolution.cxx:801
 AliTRDclusterResolution.cxx:802
 AliTRDclusterResolution.cxx:803
 AliTRDclusterResolution.cxx:804
 AliTRDclusterResolution.cxx:805
 AliTRDclusterResolution.cxx:806
 AliTRDclusterResolution.cxx:807
 AliTRDclusterResolution.cxx:808
 AliTRDclusterResolution.cxx:809
 AliTRDclusterResolution.cxx:810
 AliTRDclusterResolution.cxx:811
 AliTRDclusterResolution.cxx:812
 AliTRDclusterResolution.cxx:813
 AliTRDclusterResolution.cxx:814
 AliTRDclusterResolution.cxx:815
 AliTRDclusterResolution.cxx:816
 AliTRDclusterResolution.cxx:817
 AliTRDclusterResolution.cxx:818
 AliTRDclusterResolution.cxx:819
 AliTRDclusterResolution.cxx:820
 AliTRDclusterResolution.cxx:821
 AliTRDclusterResolution.cxx:822
 AliTRDclusterResolution.cxx:823
 AliTRDclusterResolution.cxx:824
 AliTRDclusterResolution.cxx:825
 AliTRDclusterResolution.cxx:826
 AliTRDclusterResolution.cxx:827
 AliTRDclusterResolution.cxx:828
 AliTRDclusterResolution.cxx:829
 AliTRDclusterResolution.cxx:830
 AliTRDclusterResolution.cxx:831
 AliTRDclusterResolution.cxx:832
 AliTRDclusterResolution.cxx:833
 AliTRDclusterResolution.cxx:834
 AliTRDclusterResolution.cxx:835
 AliTRDclusterResolution.cxx:836
 AliTRDclusterResolution.cxx:837
 AliTRDclusterResolution.cxx:838
 AliTRDclusterResolution.cxx:839
 AliTRDclusterResolution.cxx:840
 AliTRDclusterResolution.cxx:841
 AliTRDclusterResolution.cxx:842
 AliTRDclusterResolution.cxx:843
 AliTRDclusterResolution.cxx:844
 AliTRDclusterResolution.cxx:845
 AliTRDclusterResolution.cxx:846
 AliTRDclusterResolution.cxx:847
 AliTRDclusterResolution.cxx:848
 AliTRDclusterResolution.cxx:849
 AliTRDclusterResolution.cxx:850
 AliTRDclusterResolution.cxx:851
 AliTRDclusterResolution.cxx:852
 AliTRDclusterResolution.cxx:853
 AliTRDclusterResolution.cxx:854
 AliTRDclusterResolution.cxx:855
 AliTRDclusterResolution.cxx:856
 AliTRDclusterResolution.cxx:857
 AliTRDclusterResolution.cxx:858
 AliTRDclusterResolution.cxx:859
 AliTRDclusterResolution.cxx:860
 AliTRDclusterResolution.cxx:861
 AliTRDclusterResolution.cxx:862
 AliTRDclusterResolution.cxx:863
 AliTRDclusterResolution.cxx:864
 AliTRDclusterResolution.cxx:865
 AliTRDclusterResolution.cxx:866
 AliTRDclusterResolution.cxx:867
 AliTRDclusterResolution.cxx:868
 AliTRDclusterResolution.cxx:869
 AliTRDclusterResolution.cxx:870
 AliTRDclusterResolution.cxx:871
 AliTRDclusterResolution.cxx:872
 AliTRDclusterResolution.cxx:873
 AliTRDclusterResolution.cxx:874
 AliTRDclusterResolution.cxx:875
 AliTRDclusterResolution.cxx:876
 AliTRDclusterResolution.cxx:877
 AliTRDclusterResolution.cxx:878
 AliTRDclusterResolution.cxx:879
 AliTRDclusterResolution.cxx:880
 AliTRDclusterResolution.cxx:881
 AliTRDclusterResolution.cxx:882
 AliTRDclusterResolution.cxx:883
 AliTRDclusterResolution.cxx:884
 AliTRDclusterResolution.cxx:885
 AliTRDclusterResolution.cxx:886
 AliTRDclusterResolution.cxx:887
 AliTRDclusterResolution.cxx:888
 AliTRDclusterResolution.cxx:889
 AliTRDclusterResolution.cxx:890
 AliTRDclusterResolution.cxx:891
 AliTRDclusterResolution.cxx:892
 AliTRDclusterResolution.cxx:893
 AliTRDclusterResolution.cxx:894
 AliTRDclusterResolution.cxx:895
 AliTRDclusterResolution.cxx:896
 AliTRDclusterResolution.cxx:897
 AliTRDclusterResolution.cxx:898
 AliTRDclusterResolution.cxx:899
 AliTRDclusterResolution.cxx:900
 AliTRDclusterResolution.cxx:901
 AliTRDclusterResolution.cxx:902
 AliTRDclusterResolution.cxx:903
 AliTRDclusterResolution.cxx:904
 AliTRDclusterResolution.cxx:905
 AliTRDclusterResolution.cxx:906
 AliTRDclusterResolution.cxx:907
 AliTRDclusterResolution.cxx:908
 AliTRDclusterResolution.cxx:909
 AliTRDclusterResolution.cxx:910
 AliTRDclusterResolution.cxx:911
 AliTRDclusterResolution.cxx:912
 AliTRDclusterResolution.cxx:913
 AliTRDclusterResolution.cxx:914
 AliTRDclusterResolution.cxx:915
 AliTRDclusterResolution.cxx:916
 AliTRDclusterResolution.cxx:917
 AliTRDclusterResolution.cxx:918
 AliTRDclusterResolution.cxx:919
 AliTRDclusterResolution.cxx:920
 AliTRDclusterResolution.cxx:921
 AliTRDclusterResolution.cxx:922
 AliTRDclusterResolution.cxx:923
 AliTRDclusterResolution.cxx:924
 AliTRDclusterResolution.cxx:925
 AliTRDclusterResolution.cxx:926
 AliTRDclusterResolution.cxx:927
 AliTRDclusterResolution.cxx:928
 AliTRDclusterResolution.cxx:929
 AliTRDclusterResolution.cxx:930
 AliTRDclusterResolution.cxx:931
 AliTRDclusterResolution.cxx:932
 AliTRDclusterResolution.cxx:933
 AliTRDclusterResolution.cxx:934
 AliTRDclusterResolution.cxx:935
 AliTRDclusterResolution.cxx:936
 AliTRDclusterResolution.cxx:937
 AliTRDclusterResolution.cxx:938
 AliTRDclusterResolution.cxx:939
 AliTRDclusterResolution.cxx:940
 AliTRDclusterResolution.cxx:941
 AliTRDclusterResolution.cxx:942
 AliTRDclusterResolution.cxx:943
 AliTRDclusterResolution.cxx:944
 AliTRDclusterResolution.cxx:945
 AliTRDclusterResolution.cxx:946
 AliTRDclusterResolution.cxx:947
 AliTRDclusterResolution.cxx:948
 AliTRDclusterResolution.cxx:949
 AliTRDclusterResolution.cxx:950
 AliTRDclusterResolution.cxx:951
 AliTRDclusterResolution.cxx:952
 AliTRDclusterResolution.cxx:953
 AliTRDclusterResolution.cxx:954
 AliTRDclusterResolution.cxx:955
 AliTRDclusterResolution.cxx:956
 AliTRDclusterResolution.cxx:957
 AliTRDclusterResolution.cxx:958
 AliTRDclusterResolution.cxx:959
 AliTRDclusterResolution.cxx:960
 AliTRDclusterResolution.cxx:961
 AliTRDclusterResolution.cxx:962
 AliTRDclusterResolution.cxx:963
 AliTRDclusterResolution.cxx:964
 AliTRDclusterResolution.cxx:965
 AliTRDclusterResolution.cxx:966
 AliTRDclusterResolution.cxx:967
 AliTRDclusterResolution.cxx:968
 AliTRDclusterResolution.cxx:969
 AliTRDclusterResolution.cxx:970
 AliTRDclusterResolution.cxx:971
 AliTRDclusterResolution.cxx:972
 AliTRDclusterResolution.cxx:973
 AliTRDclusterResolution.cxx:974
 AliTRDclusterResolution.cxx:975
 AliTRDclusterResolution.cxx:976
 AliTRDclusterResolution.cxx:977
 AliTRDclusterResolution.cxx:978
 AliTRDclusterResolution.cxx:979
 AliTRDclusterResolution.cxx:980
 AliTRDclusterResolution.cxx:981
 AliTRDclusterResolution.cxx:982
 AliTRDclusterResolution.cxx:983
 AliTRDclusterResolution.cxx:984
 AliTRDclusterResolution.cxx:985
 AliTRDclusterResolution.cxx:986
 AliTRDclusterResolution.cxx:987
 AliTRDclusterResolution.cxx:988
 AliTRDclusterResolution.cxx:989
 AliTRDclusterResolution.cxx:990
 AliTRDclusterResolution.cxx:991
 AliTRDclusterResolution.cxx:992
 AliTRDclusterResolution.cxx:993
 AliTRDclusterResolution.cxx:994
 AliTRDclusterResolution.cxx:995
 AliTRDclusterResolution.cxx:996
 AliTRDclusterResolution.cxx:997
 AliTRDclusterResolution.cxx:998
 AliTRDclusterResolution.cxx:999
 AliTRDclusterResolution.cxx:1000
 AliTRDclusterResolution.cxx:1001
 AliTRDclusterResolution.cxx:1002
 AliTRDclusterResolution.cxx:1003
 AliTRDclusterResolution.cxx:1004
 AliTRDclusterResolution.cxx:1005
 AliTRDclusterResolution.cxx:1006
 AliTRDclusterResolution.cxx:1007
 AliTRDclusterResolution.cxx:1008
 AliTRDclusterResolution.cxx:1009
 AliTRDclusterResolution.cxx:1010
 AliTRDclusterResolution.cxx:1011
 AliTRDclusterResolution.cxx:1012
 AliTRDclusterResolution.cxx:1013
 AliTRDclusterResolution.cxx:1014
 AliTRDclusterResolution.cxx:1015
 AliTRDclusterResolution.cxx:1016
 AliTRDclusterResolution.cxx:1017
 AliTRDclusterResolution.cxx:1018
 AliTRDclusterResolution.cxx:1019
 AliTRDclusterResolution.cxx:1020
 AliTRDclusterResolution.cxx:1021
 AliTRDclusterResolution.cxx:1022
 AliTRDclusterResolution.cxx:1023
 AliTRDclusterResolution.cxx:1024
 AliTRDclusterResolution.cxx:1025
 AliTRDclusterResolution.cxx:1026
 AliTRDclusterResolution.cxx:1027
 AliTRDclusterResolution.cxx:1028
 AliTRDclusterResolution.cxx:1029
 AliTRDclusterResolution.cxx:1030
 AliTRDclusterResolution.cxx:1031
 AliTRDclusterResolution.cxx:1032
 AliTRDclusterResolution.cxx:1033
 AliTRDclusterResolution.cxx:1034
 AliTRDclusterResolution.cxx:1035
 AliTRDclusterResolution.cxx:1036
 AliTRDclusterResolution.cxx:1037
 AliTRDclusterResolution.cxx:1038
 AliTRDclusterResolution.cxx:1039
 AliTRDclusterResolution.cxx:1040
 AliTRDclusterResolution.cxx:1041
 AliTRDclusterResolution.cxx:1042
 AliTRDclusterResolution.cxx:1043
 AliTRDclusterResolution.cxx:1044
 AliTRDclusterResolution.cxx:1045
 AliTRDclusterResolution.cxx:1046
 AliTRDclusterResolution.cxx:1047
 AliTRDclusterResolution.cxx:1048
 AliTRDclusterResolution.cxx:1049
 AliTRDclusterResolution.cxx:1050
 AliTRDclusterResolution.cxx:1051
 AliTRDclusterResolution.cxx:1052
 AliTRDclusterResolution.cxx:1053
 AliTRDclusterResolution.cxx:1054
 AliTRDclusterResolution.cxx:1055
 AliTRDclusterResolution.cxx:1056
 AliTRDclusterResolution.cxx:1057
 AliTRDclusterResolution.cxx:1058
 AliTRDclusterResolution.cxx:1059
 AliTRDclusterResolution.cxx:1060
 AliTRDclusterResolution.cxx:1061
 AliTRDclusterResolution.cxx:1062
 AliTRDclusterResolution.cxx:1063
 AliTRDclusterResolution.cxx:1064
 AliTRDclusterResolution.cxx:1065
 AliTRDclusterResolution.cxx:1066
 AliTRDclusterResolution.cxx:1067
 AliTRDclusterResolution.cxx:1068
 AliTRDclusterResolution.cxx:1069
 AliTRDclusterResolution.cxx:1070
 AliTRDclusterResolution.cxx:1071
 AliTRDclusterResolution.cxx:1072
 AliTRDclusterResolution.cxx:1073
 AliTRDclusterResolution.cxx:1074
 AliTRDclusterResolution.cxx:1075
 AliTRDclusterResolution.cxx:1076
 AliTRDclusterResolution.cxx:1077
 AliTRDclusterResolution.cxx:1078
 AliTRDclusterResolution.cxx:1079
 AliTRDclusterResolution.cxx:1080
 AliTRDclusterResolution.cxx:1081
 AliTRDclusterResolution.cxx:1082
 AliTRDclusterResolution.cxx:1083
 AliTRDclusterResolution.cxx:1084
 AliTRDclusterResolution.cxx:1085
 AliTRDclusterResolution.cxx:1086
 AliTRDclusterResolution.cxx:1087
 AliTRDclusterResolution.cxx:1088
 AliTRDclusterResolution.cxx:1089
 AliTRDclusterResolution.cxx:1090
 AliTRDclusterResolution.cxx:1091
 AliTRDclusterResolution.cxx:1092
 AliTRDclusterResolution.cxx:1093
 AliTRDclusterResolution.cxx:1094
 AliTRDclusterResolution.cxx:1095
 AliTRDclusterResolution.cxx:1096
 AliTRDclusterResolution.cxx:1097
 AliTRDclusterResolution.cxx:1098
 AliTRDclusterResolution.cxx:1099
 AliTRDclusterResolution.cxx:1100
 AliTRDclusterResolution.cxx:1101
 AliTRDclusterResolution.cxx:1102
 AliTRDclusterResolution.cxx:1103
 AliTRDclusterResolution.cxx:1104
 AliTRDclusterResolution.cxx:1105
 AliTRDclusterResolution.cxx:1106
 AliTRDclusterResolution.cxx:1107
 AliTRDclusterResolution.cxx:1108
 AliTRDclusterResolution.cxx:1109
 AliTRDclusterResolution.cxx:1110
 AliTRDclusterResolution.cxx:1111
 AliTRDclusterResolution.cxx:1112
 AliTRDclusterResolution.cxx:1113
 AliTRDclusterResolution.cxx:1114
 AliTRDclusterResolution.cxx:1115
 AliTRDclusterResolution.cxx:1116
 AliTRDclusterResolution.cxx:1117
 AliTRDclusterResolution.cxx:1118
 AliTRDclusterResolution.cxx:1119
 AliTRDclusterResolution.cxx:1120
 AliTRDclusterResolution.cxx:1121
 AliTRDclusterResolution.cxx:1122
 AliTRDclusterResolution.cxx:1123
 AliTRDclusterResolution.cxx:1124
 AliTRDclusterResolution.cxx:1125
 AliTRDclusterResolution.cxx:1126
 AliTRDclusterResolution.cxx:1127
 AliTRDclusterResolution.cxx:1128
 AliTRDclusterResolution.cxx:1129
 AliTRDclusterResolution.cxx:1130
 AliTRDclusterResolution.cxx:1131
 AliTRDclusterResolution.cxx:1132
 AliTRDclusterResolution.cxx:1133
 AliTRDclusterResolution.cxx:1134
 AliTRDclusterResolution.cxx:1135
 AliTRDclusterResolution.cxx:1136
 AliTRDclusterResolution.cxx:1137
 AliTRDclusterResolution.cxx:1138
 AliTRDclusterResolution.cxx:1139
 AliTRDclusterResolution.cxx:1140
 AliTRDclusterResolution.cxx:1141
 AliTRDclusterResolution.cxx:1142
 AliTRDclusterResolution.cxx:1143
 AliTRDclusterResolution.cxx:1144
 AliTRDclusterResolution.cxx:1145
 AliTRDclusterResolution.cxx:1146
 AliTRDclusterResolution.cxx:1147
 AliTRDclusterResolution.cxx:1148
 AliTRDclusterResolution.cxx:1149
 AliTRDclusterResolution.cxx:1150
 AliTRDclusterResolution.cxx:1151
 AliTRDclusterResolution.cxx:1152
 AliTRDclusterResolution.cxx:1153
 AliTRDclusterResolution.cxx:1154
 AliTRDclusterResolution.cxx:1155
 AliTRDclusterResolution.cxx:1156
 AliTRDclusterResolution.cxx:1157
 AliTRDclusterResolution.cxx:1158
 AliTRDclusterResolution.cxx:1159
 AliTRDclusterResolution.cxx:1160
 AliTRDclusterResolution.cxx:1161
 AliTRDclusterResolution.cxx:1162
 AliTRDclusterResolution.cxx:1163
 AliTRDclusterResolution.cxx:1164
 AliTRDclusterResolution.cxx:1165
 AliTRDclusterResolution.cxx:1166
 AliTRDclusterResolution.cxx:1167
 AliTRDclusterResolution.cxx:1168
 AliTRDclusterResolution.cxx:1169
 AliTRDclusterResolution.cxx:1170
 AliTRDclusterResolution.cxx:1171
 AliTRDclusterResolution.cxx:1172
 AliTRDclusterResolution.cxx:1173
 AliTRDclusterResolution.cxx:1174
 AliTRDclusterResolution.cxx:1175
 AliTRDclusterResolution.cxx:1176
 AliTRDclusterResolution.cxx:1177
 AliTRDclusterResolution.cxx:1178
 AliTRDclusterResolution.cxx:1179
 AliTRDclusterResolution.cxx:1180
 AliTRDclusterResolution.cxx:1181
 AliTRDclusterResolution.cxx:1182
 AliTRDclusterResolution.cxx:1183
 AliTRDclusterResolution.cxx:1184
 AliTRDclusterResolution.cxx:1185
 AliTRDclusterResolution.cxx:1186
 AliTRDclusterResolution.cxx:1187
 AliTRDclusterResolution.cxx:1188
 AliTRDclusterResolution.cxx:1189
 AliTRDclusterResolution.cxx:1190
 AliTRDclusterResolution.cxx:1191
 AliTRDclusterResolution.cxx:1192
 AliTRDclusterResolution.cxx:1193
 AliTRDclusterResolution.cxx:1194
 AliTRDclusterResolution.cxx:1195
 AliTRDclusterResolution.cxx:1196
 AliTRDclusterResolution.cxx:1197
 AliTRDclusterResolution.cxx:1198
 AliTRDclusterResolution.cxx:1199
 AliTRDclusterResolution.cxx:1200
 AliTRDclusterResolution.cxx:1201
 AliTRDclusterResolution.cxx:1202
 AliTRDclusterResolution.cxx:1203
 AliTRDclusterResolution.cxx:1204
 AliTRDclusterResolution.cxx:1205
 AliTRDclusterResolution.cxx:1206
 AliTRDclusterResolution.cxx:1207
 AliTRDclusterResolution.cxx:1208
 AliTRDclusterResolution.cxx:1209
 AliTRDclusterResolution.cxx:1210
 AliTRDclusterResolution.cxx:1211
 AliTRDclusterResolution.cxx:1212
 AliTRDclusterResolution.cxx:1213
 AliTRDclusterResolution.cxx:1214
 AliTRDclusterResolution.cxx:1215
 AliTRDclusterResolution.cxx:1216
 AliTRDclusterResolution.cxx:1217
 AliTRDclusterResolution.cxx:1218
 AliTRDclusterResolution.cxx:1219
 AliTRDclusterResolution.cxx:1220
 AliTRDclusterResolution.cxx:1221
 AliTRDclusterResolution.cxx:1222
 AliTRDclusterResolution.cxx:1223
 AliTRDclusterResolution.cxx:1224
 AliTRDclusterResolution.cxx:1225
 AliTRDclusterResolution.cxx:1226
 AliTRDclusterResolution.cxx:1227
 AliTRDclusterResolution.cxx:1228
 AliTRDclusterResolution.cxx:1229
 AliTRDclusterResolution.cxx:1230
 AliTRDclusterResolution.cxx:1231
 AliTRDclusterResolution.cxx:1232
 AliTRDclusterResolution.cxx:1233
 AliTRDclusterResolution.cxx:1234
 AliTRDclusterResolution.cxx:1235
 AliTRDclusterResolution.cxx:1236
 AliTRDclusterResolution.cxx:1237
 AliTRDclusterResolution.cxx:1238
 AliTRDclusterResolution.cxx:1239
 AliTRDclusterResolution.cxx:1240
 AliTRDclusterResolution.cxx:1241
 AliTRDclusterResolution.cxx:1242
 AliTRDclusterResolution.cxx:1243
 AliTRDclusterResolution.cxx:1244
 AliTRDclusterResolution.cxx:1245
 AliTRDclusterResolution.cxx:1246
 AliTRDclusterResolution.cxx:1247
 AliTRDclusterResolution.cxx:1248
 AliTRDclusterResolution.cxx:1249
 AliTRDclusterResolution.cxx:1250
 AliTRDclusterResolution.cxx:1251
 AliTRDclusterResolution.cxx:1252
 AliTRDclusterResolution.cxx:1253
 AliTRDclusterResolution.cxx:1254
 AliTRDclusterResolution.cxx:1255
 AliTRDclusterResolution.cxx:1256
 AliTRDclusterResolution.cxx:1257
 AliTRDclusterResolution.cxx:1258
 AliTRDclusterResolution.cxx:1259
 AliTRDclusterResolution.cxx:1260
 AliTRDclusterResolution.cxx:1261
 AliTRDclusterResolution.cxx:1262
 AliTRDclusterResolution.cxx:1263
 AliTRDclusterResolution.cxx:1264
 AliTRDclusterResolution.cxx:1265
 AliTRDclusterResolution.cxx:1266
 AliTRDclusterResolution.cxx:1267
 AliTRDclusterResolution.cxx:1268
 AliTRDclusterResolution.cxx:1269
 AliTRDclusterResolution.cxx:1270
 AliTRDclusterResolution.cxx:1271
 AliTRDclusterResolution.cxx:1272
 AliTRDclusterResolution.cxx:1273
 AliTRDclusterResolution.cxx:1274
 AliTRDclusterResolution.cxx:1275
 AliTRDclusterResolution.cxx:1276
 AliTRDclusterResolution.cxx:1277
 AliTRDclusterResolution.cxx:1278
 AliTRDclusterResolution.cxx:1279
 AliTRDclusterResolution.cxx:1280
 AliTRDclusterResolution.cxx:1281
 AliTRDclusterResolution.cxx:1282
 AliTRDclusterResolution.cxx:1283
 AliTRDclusterResolution.cxx:1284
 AliTRDclusterResolution.cxx:1285
 AliTRDclusterResolution.cxx:1286
 AliTRDclusterResolution.cxx:1287
 AliTRDclusterResolution.cxx:1288
 AliTRDclusterResolution.cxx:1289
 AliTRDclusterResolution.cxx:1290
 AliTRDclusterResolution.cxx:1291
 AliTRDclusterResolution.cxx:1292
 AliTRDclusterResolution.cxx:1293
 AliTRDclusterResolution.cxx:1294
 AliTRDclusterResolution.cxx:1295
 AliTRDclusterResolution.cxx:1296
 AliTRDclusterResolution.cxx:1297
 AliTRDclusterResolution.cxx:1298
 AliTRDclusterResolution.cxx:1299
 AliTRDclusterResolution.cxx:1300
 AliTRDclusterResolution.cxx:1301
 AliTRDclusterResolution.cxx:1302
 AliTRDclusterResolution.cxx:1303
 AliTRDclusterResolution.cxx:1304
 AliTRDclusterResolution.cxx:1305
 AliTRDclusterResolution.cxx:1306
 AliTRDclusterResolution.cxx:1307
 AliTRDclusterResolution.cxx:1308
 AliTRDclusterResolution.cxx:1309
 AliTRDclusterResolution.cxx:1310
 AliTRDclusterResolution.cxx:1311
 AliTRDclusterResolution.cxx:1312
 AliTRDclusterResolution.cxx:1313
 AliTRDclusterResolution.cxx:1314
 AliTRDclusterResolution.cxx:1315
 AliTRDclusterResolution.cxx:1316
 AliTRDclusterResolution.cxx:1317
 AliTRDclusterResolution.cxx:1318
 AliTRDclusterResolution.cxx:1319
 AliTRDclusterResolution.cxx:1320
 AliTRDclusterResolution.cxx:1321
 AliTRDclusterResolution.cxx:1322
 AliTRDclusterResolution.cxx:1323
 AliTRDclusterResolution.cxx:1324
 AliTRDclusterResolution.cxx:1325
 AliTRDclusterResolution.cxx:1326
 AliTRDclusterResolution.cxx:1327
 AliTRDclusterResolution.cxx:1328
 AliTRDclusterResolution.cxx:1329
 AliTRDclusterResolution.cxx:1330
 AliTRDclusterResolution.cxx:1331
 AliTRDclusterResolution.cxx:1332
 AliTRDclusterResolution.cxx:1333
 AliTRDclusterResolution.cxx:1334
 AliTRDclusterResolution.cxx:1335
 AliTRDclusterResolution.cxx:1336
 AliTRDclusterResolution.cxx:1337
 AliTRDclusterResolution.cxx:1338
 AliTRDclusterResolution.cxx:1339
 AliTRDclusterResolution.cxx:1340
 AliTRDclusterResolution.cxx:1341
 AliTRDclusterResolution.cxx:1342
 AliTRDclusterResolution.cxx:1343
 AliTRDclusterResolution.cxx:1344
 AliTRDclusterResolution.cxx:1345
 AliTRDclusterResolution.cxx:1346
 AliTRDclusterResolution.cxx:1347
 AliTRDclusterResolution.cxx:1348
 AliTRDclusterResolution.cxx:1349
 AliTRDclusterResolution.cxx:1350
 AliTRDclusterResolution.cxx:1351
 AliTRDclusterResolution.cxx:1352
 AliTRDclusterResolution.cxx:1353
 AliTRDclusterResolution.cxx:1354
 AliTRDclusterResolution.cxx:1355
 AliTRDclusterResolution.cxx:1356
 AliTRDclusterResolution.cxx:1357
 AliTRDclusterResolution.cxx:1358
 AliTRDclusterResolution.cxx:1359
 AliTRDclusterResolution.cxx:1360
 AliTRDclusterResolution.cxx:1361
 AliTRDclusterResolution.cxx:1362
 AliTRDclusterResolution.cxx:1363
 AliTRDclusterResolution.cxx:1364
 AliTRDclusterResolution.cxx:1365
 AliTRDclusterResolution.cxx:1366
 AliTRDclusterResolution.cxx:1367
 AliTRDclusterResolution.cxx:1368
 AliTRDclusterResolution.cxx:1369
 AliTRDclusterResolution.cxx:1370
 AliTRDclusterResolution.cxx:1371
 AliTRDclusterResolution.cxx:1372
 AliTRDclusterResolution.cxx:1373
 AliTRDclusterResolution.cxx:1374
 AliTRDclusterResolution.cxx:1375
 AliTRDclusterResolution.cxx:1376
 AliTRDclusterResolution.cxx:1377
 AliTRDclusterResolution.cxx:1378
 AliTRDclusterResolution.cxx:1379
 AliTRDclusterResolution.cxx:1380
 AliTRDclusterResolution.cxx:1381
 AliTRDclusterResolution.cxx:1382
 AliTRDclusterResolution.cxx:1383
 AliTRDclusterResolution.cxx:1384
 AliTRDclusterResolution.cxx:1385
 AliTRDclusterResolution.cxx:1386
 AliTRDclusterResolution.cxx:1387
 AliTRDclusterResolution.cxx:1388
 AliTRDclusterResolution.cxx:1389
 AliTRDclusterResolution.cxx:1390
 AliTRDclusterResolution.cxx:1391
 AliTRDclusterResolution.cxx:1392
 AliTRDclusterResolution.cxx:1393
 AliTRDclusterResolution.cxx:1394
 AliTRDclusterResolution.cxx:1395
 AliTRDclusterResolution.cxx:1396
 AliTRDclusterResolution.cxx:1397
 AliTRDclusterResolution.cxx:1398
 AliTRDclusterResolution.cxx:1399
 AliTRDclusterResolution.cxx:1400
 AliTRDclusterResolution.cxx:1401
 AliTRDclusterResolution.cxx:1402
 AliTRDclusterResolution.cxx:1403
 AliTRDclusterResolution.cxx:1404
 AliTRDclusterResolution.cxx:1405
 AliTRDclusterResolution.cxx:1406
 AliTRDclusterResolution.cxx:1407
 AliTRDclusterResolution.cxx:1408
 AliTRDclusterResolution.cxx:1409
 AliTRDclusterResolution.cxx:1410
 AliTRDclusterResolution.cxx:1411
 AliTRDclusterResolution.cxx:1412
 AliTRDclusterResolution.cxx:1413
 AliTRDclusterResolution.cxx:1414
 AliTRDclusterResolution.cxx:1415
 AliTRDclusterResolution.cxx:1416
 AliTRDclusterResolution.cxx:1417
 AliTRDclusterResolution.cxx:1418
 AliTRDclusterResolution.cxx:1419
 AliTRDclusterResolution.cxx:1420
 AliTRDclusterResolution.cxx:1421
 AliTRDclusterResolution.cxx:1422
 AliTRDclusterResolution.cxx:1423
 AliTRDclusterResolution.cxx:1424
 AliTRDclusterResolution.cxx:1425
 AliTRDclusterResolution.cxx:1426
 AliTRDclusterResolution.cxx:1427
 AliTRDclusterResolution.cxx:1428
 AliTRDclusterResolution.cxx:1429
 AliTRDclusterResolution.cxx:1430
 AliTRDclusterResolution.cxx:1431
 AliTRDclusterResolution.cxx:1432
 AliTRDclusterResolution.cxx:1433
 AliTRDclusterResolution.cxx:1434
 AliTRDclusterResolution.cxx:1435
 AliTRDclusterResolution.cxx:1436
 AliTRDclusterResolution.cxx:1437
 AliTRDclusterResolution.cxx:1438
 AliTRDclusterResolution.cxx:1439
 AliTRDclusterResolution.cxx:1440
 AliTRDclusterResolution.cxx:1441
 AliTRDclusterResolution.cxx:1442
 AliTRDclusterResolution.cxx:1443
 AliTRDclusterResolution.cxx:1444
 AliTRDclusterResolution.cxx:1445
 AliTRDclusterResolution.cxx:1446
 AliTRDclusterResolution.cxx:1447
 AliTRDclusterResolution.cxx:1448
 AliTRDclusterResolution.cxx:1449
 AliTRDclusterResolution.cxx:1450
 AliTRDclusterResolution.cxx:1451
 AliTRDclusterResolution.cxx:1452
 AliTRDclusterResolution.cxx:1453
 AliTRDclusterResolution.cxx:1454
 AliTRDclusterResolution.cxx:1455
 AliTRDclusterResolution.cxx:1456
 AliTRDclusterResolution.cxx:1457
 AliTRDclusterResolution.cxx:1458
 AliTRDclusterResolution.cxx:1459
 AliTRDclusterResolution.cxx:1460
 AliTRDclusterResolution.cxx:1461
 AliTRDclusterResolution.cxx:1462
 AliTRDclusterResolution.cxx:1463
 AliTRDclusterResolution.cxx:1464
 AliTRDclusterResolution.cxx:1465
 AliTRDclusterResolution.cxx:1466
 AliTRDclusterResolution.cxx:1467
 AliTRDclusterResolution.cxx:1468
 AliTRDclusterResolution.cxx:1469
 AliTRDclusterResolution.cxx:1470
 AliTRDclusterResolution.cxx:1471
 AliTRDclusterResolution.cxx:1472
 AliTRDclusterResolution.cxx:1473
 AliTRDclusterResolution.cxx:1474
 AliTRDclusterResolution.cxx:1475
 AliTRDclusterResolution.cxx:1476
 AliTRDclusterResolution.cxx:1477
 AliTRDclusterResolution.cxx:1478
 AliTRDclusterResolution.cxx:1479
 AliTRDclusterResolution.cxx:1480
 AliTRDclusterResolution.cxx:1481
 AliTRDclusterResolution.cxx:1482
 AliTRDclusterResolution.cxx:1483
 AliTRDclusterResolution.cxx:1484
 AliTRDclusterResolution.cxx:1485
 AliTRDclusterResolution.cxx:1486
 AliTRDclusterResolution.cxx:1487
 AliTRDclusterResolution.cxx:1488
 AliTRDclusterResolution.cxx:1489
 AliTRDclusterResolution.cxx:1490
 AliTRDclusterResolution.cxx:1491
 AliTRDclusterResolution.cxx:1492
 AliTRDclusterResolution.cxx:1493
 AliTRDclusterResolution.cxx:1494
 AliTRDclusterResolution.cxx:1495
 AliTRDclusterResolution.cxx:1496
 AliTRDclusterResolution.cxx:1497
 AliTRDclusterResolution.cxx:1498
 AliTRDclusterResolution.cxx:1499
 AliTRDclusterResolution.cxx:1500
 AliTRDclusterResolution.cxx:1501
 AliTRDclusterResolution.cxx:1502
 AliTRDclusterResolution.cxx:1503
 AliTRDclusterResolution.cxx:1504
 AliTRDclusterResolution.cxx:1505
 AliTRDclusterResolution.cxx:1506
 AliTRDclusterResolution.cxx:1507
 AliTRDclusterResolution.cxx:1508
 AliTRDclusterResolution.cxx:1509
 AliTRDclusterResolution.cxx:1510
 AliTRDclusterResolution.cxx:1511