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-commercial 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$ */


//-----------------------------------------------------------------------------
//
//
//
//  Origin:   Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
//
//  Declaration of class AliTPCRF1D
//
//-----------------------------------------------------------------------------

//

#include <RVersion.h>
#include <Riostream.h>
#include <TCanvas.h>
#include <TClass.h>
#include <TF2.h>
#include <TH1.h>
#include <TMath.h>
#include <TPad.h>
#include <TString.h>
#include <TStyle.h>

#include "AliTPCRF1D.h"

extern TStyle * gStyle; 

Int_t   AliTPCRF1D::fgNRF=100;  //default  number of interpolation points
Float_t AliTPCRF1D::fgRFDSTEP=0.01; //default step in cm

static Double_t funGauss(Double_t *x, Double_t * par)
{
  //Gauss function  -needde by the generic function object 
  return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
}

static Double_t funCosh(Double_t *x, Double_t * par)
{
  //Cosh function  -needde by the generic function object 
  return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));  
}    

static Double_t funGati(Double_t *x, Double_t * par)
{
  //Gati function  -needde by the generic function object 
  Float_t k3=par[1];
  Float_t k3R=TMath::Sqrt(k3);
  Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
  Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
  Float_t l=x[0]/par[0];
  Float_t tan2=TMath::TanH(k2*l);
  tan2*=tan2;
  Float_t res = k1*(1-tan2)/(1+k3*tan2);  
  return res;  
}    

///////////////////////////////////////////////////////////////////////////
///////////////////////////////////////////////////////////////////////////

ClassImp(AliTPCRF1D)


AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
           :TObject(),
	    fNRF(0),
            fDSTEPM1(0.),
            fcharge(0),
            forigsigma(0.),
            fpadWidth(3.5),
            fkNorm(0.5),
            fInteg(0.),
            fGRF(0),
            fSigma(0.),
            fOffset(0.),
            fDirect(kFALSE),
            fPadDistance(0.)
{
  //default constructor for response function object
  fDirect=direct;
  if (np!=0)fNRF = np;
  else (fNRF=fgNRF);
  fcharge = new Float_t[fNRF];
  if (step>0) fDSTEPM1=1./step;
  else fDSTEPM1 = 1./fgRFDSTEP;
  for(Int_t i=0;i<5;i++) {
    funParam[i]=0.;
    fType[i]=0;
  }
  
}

AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
           :TObject(prf),
	    fNRF(prf.fNRF),
            fDSTEPM1(prf.fDSTEPM1),
            fcharge(0),
            forigsigma(prf.forigsigma),
            fpadWidth(prf.fpadWidth),
            fkNorm(prf.fkNorm),
            fInteg(prf.fInteg),
            fGRF(new TF1(*(prf.fGRF))),
            fSigma(prf.fSigma),
            fOffset(prf.fOffset),
            fDirect(prf.fDirect),
            fPadDistance(prf.fPadDistance)
{
  //
  //
  for(Int_t i=0;i<5;i++) {
    funParam[i]=0.;
    fType[i]=0;
  }
  fcharge = new Float_t[fNRF];
  memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));

  //PH Change the name (add 0 to the end)
  TString s(fGRF->GetName());
  s+="0";
  fGRF->SetName(s.Data());
}

AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
{
  if(this!=&prf) {
    TObject::operator=(prf);
    fNRF=prf.fNRF;
    fDSTEPM1=prf.fDSTEPM1;
    delete [] fcharge;
    fcharge = new Float_t[fNRF];
    memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
    forigsigma=prf.forigsigma;
    fpadWidth=prf.fpadWidth;
    fkNorm=prf.fkNorm;
    fInteg=prf.fInteg;
    delete fGRF;
    fGRF=new TF1(*(prf.fGRF));
   //PH Change the name (add 0 to the end)
    TString s(fGRF->GetName());
    s+="0";
    fGRF->SetName(s.Data());
    fSigma=prf.fSigma;
    fOffset=prf.fOffset;
    fDirect=prf.fDirect;
    fPadDistance=prf.fPadDistance;
  }
  return *this;
}



AliTPCRF1D::~AliTPCRF1D()
{
  //
  delete [] fcharge;
  delete fGRF;
}

Float_t AliTPCRF1D::GetRF(Float_t xin)
{
  //function which return response
  //for the charge in distance xin 
  //return linear aproximation of RF
  Float_t x = (xin-fOffset)*fDSTEPM1+fNRF/2;
  Int_t i1=Int_t(x);
  if (x<0) i1-=1;
  Float_t res=0;
  if (i1+1<fNRF &&i1>0)
    res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));    
  return res;
}

Float_t  AliTPCRF1D::GetGRF(Float_t xin)
{  
  //function which returnoriginal charge distribution
  //this function is just normalised for fKnorm
  if (fGRF != 0 ) 
    return fkNorm*fGRF->Eval(xin)/fInteg;
      else
    return 0.;
}

   
void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
		       Float_t kNorm, Float_t sigma)
{
  //adjust parameters of the original charge distribution
  //and pad size parameters
   fpadWidth = padwidth;
   fGRF = GRF;
   fkNorm = kNorm;
   if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
   forigsigma=sigma;
   fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12); 
   //sprintf(fType,"User");
   snprintf(fType,5,"User");
   //   Update();   
}
  

void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
		      Float_t kNorm)
{
  // 
  // set parameters for Gauss generic charge distribution
  //
  fpadWidth = padWidth;
  fkNorm = kNorm;
  if (fGRF !=0 ) fGRF->Delete();
  fGRF = new TF1("funGauss",funGauss,-5,5,1);
  funParam[0]=sigma;
  forigsigma=sigma;
  fGRF->SetParameters(funParam);
   fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12); 
  //by default I set the step as one tenth of sigma  
  //sprintf(fType,"Gauss");
   snprintf(fType,5,"Gauss");
}

void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
		     Float_t kNorm)
{
  // 
  // set parameters for Cosh generic charge distribution
  //
  fpadWidth = padWidth;
  fkNorm = kNorm;
  if (fGRF !=0 ) fGRF->Delete();
  fGRF = new TF1("funCosh",	funCosh, -5.,5.,2);   
  funParam[0]=sigma;
  fGRF->SetParameters(funParam);
  forigsigma=sigma;
  fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12); 
  //by default I set the step as one tenth of sigma
  //sprintf(fType,"Cosh");
  snprintf(fType,5,"Cosh");
}

void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
		     Float_t kNorm)
{
  // 
  // set parameters for Gati generic charge distribution
  //
  fpadWidth = padWidth;
  fkNorm = kNorm;
  if (fGRF !=0 ) fGRF->Delete();
  fGRF = new TF1("funGati",	funGati, -5.,5.,2);   
  funParam[0]=padDistance;
  funParam[1]=K3;  
  fGRF->SetParameters(funParam);
  forigsigma=padDistance;
  fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12); 
  //by default I set the step as one tenth of sigma
  //sprintf(fType,"Gati");
  snprintf(fType,5,"Gati");
}



void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
{ 
  //
  //Draw prf in selected region <x1,x2> with nuber of diviision = n
  //
  char s[100];
  TCanvas  * c1 = new TCanvas("canRF","Pad response function",700,900);
  c1->cd();
  TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
  pad1->Draw();
  TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
  pad2->Draw();

  //sprintf(s,"RF response function for %1.2f cm pad width",
  //	  fpadWidth); 
  snprintf(s,60,"RF response function for %1.2f cm pad width",fpadWidth); 
  pad1->cd();
  TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
  pad2->cd();
   gStyle->SetOptFit(1);
   gStyle->SetOptStat(0); 
  TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
  Float_t x=x1;
  Float_t y1;
  Float_t y2;

  for (Float_t i = 0;i<N+1;i++)
    {
      x+=(x2-x1)/Float_t(N);
      y1 = GetRF(x);
      hRFc->Fill(x,y1);
      y2 = GetGRF(x);
      hRFo->Fill(x,y2);      
    };
  pad1->cd();
  hRFo->Fit("gaus");
  pad2->cd();
  hRFc->Fit("gaus");
}

void AliTPCRF1D::Update()
{
  //
  //update fields  with interpolated values for
  //PRF calculation

  //at the begining initialize to 0
  for (Int_t i =0; i<fNRF;i++)  fcharge[i] = 0;
  if ( fGRF == 0 ) return;
  // This form is no longer available 
#if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
  fInteg  = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
#else
  TArrayD savParam(fGRF->GetNpar(), fGRF->GetParameters());
  fGRF->SetParameters(funParam);
  fInteg  = fGRF->Integral(-5*forigsigma,5*forigsigma,0.00001);
#endif
  if ( fInteg == 0 ) fInteg = 1; 
  if (fDirect==kFALSE){
  //integrate charge over pad for different distance of pad
  for (Int_t i =0; i<fNRF;i++)
    {      //x in cm fpadWidth in cm
      Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
      Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
      Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
#if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
      fcharge[i] = fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
#else
      fcharge[i] = fkNorm*fGRF->Integral(x1,x2,0.0001)/fInteg;
#endif
    };   
  }
  else{
    for (Int_t i =0; i<fNRF;i++)
      {      //x in cm fpadWidth in cm
	Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
	fcharge[i] = fkNorm*fGRF->Eval(x);
      };   
  }  
  fSigma = 0; 
  Float_t sum =0;
  Float_t mean=0;
  for (Float_t  x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
    {      //x in cm fpadWidth in cm
      Float_t weight = GetRF(x+fOffset);
      fSigma+=x*x*weight; 
      mean+=x*weight;
      sum+=weight;
    };  
  if (sum>0){
    mean/=sum;
    fSigma = TMath::Sqrt(fSigma/sum-mean*mean);   
  }
  else fSigma=0; 
#if ROOT_VERSION_CODE >= ROOT_VERSION(5,99,0)
  fGRF->SetParameters(savParam.GetArray());
#endif
}

void AliTPCRF1D::Streamer(TBuffer &R__b)
{
   // Stream an object of class AliTPCRF1D.
   if (R__b.IsReading()) {
      AliTPCRF1D::Class()->ReadBuffer(R__b, this);
      //read functions
 
      if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("funGauss",funGauss,-5.,5.,4);}
      if (strncmp(fType,"Cosh",3)==0)  {delete fGRF; fGRF = new TF1("funCosh",funCosh,-5.,5.,4);}
      if (strncmp(fType,"Gati",3)==0)  {delete fGRF; fGRF = new TF1("funGati",funGati,-5.,5.,4);}  
      if (fGRF) fGRF->SetParameters(funParam);     

   } else {
      AliTPCRF1D::Class()->WriteBuffer(R__b, this);
   }
}


Double_t  AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
  //
  // Gamma 4 Time response function of ALTRO
  //
  if (x<0) return 0;
  Double_t g1 = TMath::Exp(-4.*x/p1);
  Double_t g2 = TMath::Power(x/p1,4);
  return p0*g1*g2;
}
 
 AliTPCRF1D.cxx:1
 AliTPCRF1D.cxx:2
 AliTPCRF1D.cxx:3
 AliTPCRF1D.cxx:4
 AliTPCRF1D.cxx:5
 AliTPCRF1D.cxx:6
 AliTPCRF1D.cxx:7
 AliTPCRF1D.cxx:8
 AliTPCRF1D.cxx:9
 AliTPCRF1D.cxx:10
 AliTPCRF1D.cxx:11
 AliTPCRF1D.cxx:12
 AliTPCRF1D.cxx:13
 AliTPCRF1D.cxx:14
 AliTPCRF1D.cxx:15
 AliTPCRF1D.cxx:16
 AliTPCRF1D.cxx:17
 AliTPCRF1D.cxx:18
 AliTPCRF1D.cxx:19
 AliTPCRF1D.cxx:20
 AliTPCRF1D.cxx:21
 AliTPCRF1D.cxx:22
 AliTPCRF1D.cxx:23
 AliTPCRF1D.cxx:24
 AliTPCRF1D.cxx:25
 AliTPCRF1D.cxx:26
 AliTPCRF1D.cxx:27
 AliTPCRF1D.cxx:28
 AliTPCRF1D.cxx:29
 AliTPCRF1D.cxx:30
 AliTPCRF1D.cxx:31
 AliTPCRF1D.cxx:32
 AliTPCRF1D.cxx:33
 AliTPCRF1D.cxx:34
 AliTPCRF1D.cxx:35
 AliTPCRF1D.cxx:36
 AliTPCRF1D.cxx:37
 AliTPCRF1D.cxx:38
 AliTPCRF1D.cxx:39
 AliTPCRF1D.cxx:40
 AliTPCRF1D.cxx:41
 AliTPCRF1D.cxx:42
 AliTPCRF1D.cxx:43
 AliTPCRF1D.cxx:44
 AliTPCRF1D.cxx:45
 AliTPCRF1D.cxx:46
 AliTPCRF1D.cxx:47
 AliTPCRF1D.cxx:48
 AliTPCRF1D.cxx:49
 AliTPCRF1D.cxx:50
 AliTPCRF1D.cxx:51
 AliTPCRF1D.cxx:52
 AliTPCRF1D.cxx:53
 AliTPCRF1D.cxx:54
 AliTPCRF1D.cxx:55
 AliTPCRF1D.cxx:56
 AliTPCRF1D.cxx:57
 AliTPCRF1D.cxx:58
 AliTPCRF1D.cxx:59
 AliTPCRF1D.cxx:60
 AliTPCRF1D.cxx:61
 AliTPCRF1D.cxx:62
 AliTPCRF1D.cxx:63
 AliTPCRF1D.cxx:64
 AliTPCRF1D.cxx:65
 AliTPCRF1D.cxx:66
 AliTPCRF1D.cxx:67
 AliTPCRF1D.cxx:68
 AliTPCRF1D.cxx:69
 AliTPCRF1D.cxx:70
 AliTPCRF1D.cxx:71
 AliTPCRF1D.cxx:72
 AliTPCRF1D.cxx:73
 AliTPCRF1D.cxx:74
 AliTPCRF1D.cxx:75
 AliTPCRF1D.cxx:76
 AliTPCRF1D.cxx:77
 AliTPCRF1D.cxx:78
 AliTPCRF1D.cxx:79
 AliTPCRF1D.cxx:80
 AliTPCRF1D.cxx:81
 AliTPCRF1D.cxx:82
 AliTPCRF1D.cxx:83
 AliTPCRF1D.cxx:84
 AliTPCRF1D.cxx:85
 AliTPCRF1D.cxx:86
 AliTPCRF1D.cxx:87
 AliTPCRF1D.cxx:88
 AliTPCRF1D.cxx:89
 AliTPCRF1D.cxx:90
 AliTPCRF1D.cxx:91
 AliTPCRF1D.cxx:92
 AliTPCRF1D.cxx:93
 AliTPCRF1D.cxx:94
 AliTPCRF1D.cxx:95
 AliTPCRF1D.cxx:96
 AliTPCRF1D.cxx:97
 AliTPCRF1D.cxx:98
 AliTPCRF1D.cxx:99
 AliTPCRF1D.cxx:100
 AliTPCRF1D.cxx:101
 AliTPCRF1D.cxx:102
 AliTPCRF1D.cxx:103
 AliTPCRF1D.cxx:104
 AliTPCRF1D.cxx:105
 AliTPCRF1D.cxx:106
 AliTPCRF1D.cxx:107
 AliTPCRF1D.cxx:108
 AliTPCRF1D.cxx:109
 AliTPCRF1D.cxx:110
 AliTPCRF1D.cxx:111
 AliTPCRF1D.cxx:112
 AliTPCRF1D.cxx:113
 AliTPCRF1D.cxx:114
 AliTPCRF1D.cxx:115
 AliTPCRF1D.cxx:116
 AliTPCRF1D.cxx:117
 AliTPCRF1D.cxx:118
 AliTPCRF1D.cxx:119
 AliTPCRF1D.cxx:120
 AliTPCRF1D.cxx:121
 AliTPCRF1D.cxx:122
 AliTPCRF1D.cxx:123
 AliTPCRF1D.cxx:124
 AliTPCRF1D.cxx:125
 AliTPCRF1D.cxx:126
 AliTPCRF1D.cxx:127
 AliTPCRF1D.cxx:128
 AliTPCRF1D.cxx:129
 AliTPCRF1D.cxx:130
 AliTPCRF1D.cxx:131
 AliTPCRF1D.cxx:132
 AliTPCRF1D.cxx:133
 AliTPCRF1D.cxx:134
 AliTPCRF1D.cxx:135
 AliTPCRF1D.cxx:136
 AliTPCRF1D.cxx:137
 AliTPCRF1D.cxx:138
 AliTPCRF1D.cxx:139
 AliTPCRF1D.cxx:140
 AliTPCRF1D.cxx:141
 AliTPCRF1D.cxx:142
 AliTPCRF1D.cxx:143
 AliTPCRF1D.cxx:144
 AliTPCRF1D.cxx:145
 AliTPCRF1D.cxx:146
 AliTPCRF1D.cxx:147
 AliTPCRF1D.cxx:148
 AliTPCRF1D.cxx:149
 AliTPCRF1D.cxx:150
 AliTPCRF1D.cxx:151
 AliTPCRF1D.cxx:152
 AliTPCRF1D.cxx:153
 AliTPCRF1D.cxx:154
 AliTPCRF1D.cxx:155
 AliTPCRF1D.cxx:156
 AliTPCRF1D.cxx:157
 AliTPCRF1D.cxx:158
 AliTPCRF1D.cxx:159
 AliTPCRF1D.cxx:160
 AliTPCRF1D.cxx:161
 AliTPCRF1D.cxx:162
 AliTPCRF1D.cxx:163
 AliTPCRF1D.cxx:164
 AliTPCRF1D.cxx:165
 AliTPCRF1D.cxx:166
 AliTPCRF1D.cxx:167
 AliTPCRF1D.cxx:168
 AliTPCRF1D.cxx:169
 AliTPCRF1D.cxx:170
 AliTPCRF1D.cxx:171
 AliTPCRF1D.cxx:172
 AliTPCRF1D.cxx:173
 AliTPCRF1D.cxx:174
 AliTPCRF1D.cxx:175
 AliTPCRF1D.cxx:176
 AliTPCRF1D.cxx:177
 AliTPCRF1D.cxx:178
 AliTPCRF1D.cxx:179
 AliTPCRF1D.cxx:180
 AliTPCRF1D.cxx:181
 AliTPCRF1D.cxx:182
 AliTPCRF1D.cxx:183
 AliTPCRF1D.cxx:184
 AliTPCRF1D.cxx:185
 AliTPCRF1D.cxx:186
 AliTPCRF1D.cxx:187
 AliTPCRF1D.cxx:188
 AliTPCRF1D.cxx:189
 AliTPCRF1D.cxx:190
 AliTPCRF1D.cxx:191
 AliTPCRF1D.cxx:192
 AliTPCRF1D.cxx:193
 AliTPCRF1D.cxx:194
 AliTPCRF1D.cxx:195
 AliTPCRF1D.cxx:196
 AliTPCRF1D.cxx:197
 AliTPCRF1D.cxx:198
 AliTPCRF1D.cxx:199
 AliTPCRF1D.cxx:200
 AliTPCRF1D.cxx:201
 AliTPCRF1D.cxx:202
 AliTPCRF1D.cxx:203
 AliTPCRF1D.cxx:204
 AliTPCRF1D.cxx:205
 AliTPCRF1D.cxx:206
 AliTPCRF1D.cxx:207
 AliTPCRF1D.cxx:208
 AliTPCRF1D.cxx:209
 AliTPCRF1D.cxx:210
 AliTPCRF1D.cxx:211
 AliTPCRF1D.cxx:212
 AliTPCRF1D.cxx:213
 AliTPCRF1D.cxx:214
 AliTPCRF1D.cxx:215
 AliTPCRF1D.cxx:216
 AliTPCRF1D.cxx:217
 AliTPCRF1D.cxx:218
 AliTPCRF1D.cxx:219
 AliTPCRF1D.cxx:220
 AliTPCRF1D.cxx:221
 AliTPCRF1D.cxx:222
 AliTPCRF1D.cxx:223
 AliTPCRF1D.cxx:224
 AliTPCRF1D.cxx:225
 AliTPCRF1D.cxx:226
 AliTPCRF1D.cxx:227
 AliTPCRF1D.cxx:228
 AliTPCRF1D.cxx:229
 AliTPCRF1D.cxx:230
 AliTPCRF1D.cxx:231
 AliTPCRF1D.cxx:232
 AliTPCRF1D.cxx:233
 AliTPCRF1D.cxx:234
 AliTPCRF1D.cxx:235
 AliTPCRF1D.cxx:236
 AliTPCRF1D.cxx:237
 AliTPCRF1D.cxx:238
 AliTPCRF1D.cxx:239
 AliTPCRF1D.cxx:240
 AliTPCRF1D.cxx:241
 AliTPCRF1D.cxx:242
 AliTPCRF1D.cxx:243
 AliTPCRF1D.cxx:244
 AliTPCRF1D.cxx:245
 AliTPCRF1D.cxx:246
 AliTPCRF1D.cxx:247
 AliTPCRF1D.cxx:248
 AliTPCRF1D.cxx:249
 AliTPCRF1D.cxx:250
 AliTPCRF1D.cxx:251
 AliTPCRF1D.cxx:252
 AliTPCRF1D.cxx:253
 AliTPCRF1D.cxx:254
 AliTPCRF1D.cxx:255
 AliTPCRF1D.cxx:256
 AliTPCRF1D.cxx:257
 AliTPCRF1D.cxx:258
 AliTPCRF1D.cxx:259
 AliTPCRF1D.cxx:260
 AliTPCRF1D.cxx:261
 AliTPCRF1D.cxx:262
 AliTPCRF1D.cxx:263
 AliTPCRF1D.cxx:264
 AliTPCRF1D.cxx:265
 AliTPCRF1D.cxx:266
 AliTPCRF1D.cxx:267
 AliTPCRF1D.cxx:268
 AliTPCRF1D.cxx:269
 AliTPCRF1D.cxx:270
 AliTPCRF1D.cxx:271
 AliTPCRF1D.cxx:272
 AliTPCRF1D.cxx:273
 AliTPCRF1D.cxx:274
 AliTPCRF1D.cxx:275
 AliTPCRF1D.cxx:276
 AliTPCRF1D.cxx:277
 AliTPCRF1D.cxx:278
 AliTPCRF1D.cxx:279
 AliTPCRF1D.cxx:280
 AliTPCRF1D.cxx:281
 AliTPCRF1D.cxx:282
 AliTPCRF1D.cxx:283
 AliTPCRF1D.cxx:284
 AliTPCRF1D.cxx:285
 AliTPCRF1D.cxx:286
 AliTPCRF1D.cxx:287
 AliTPCRF1D.cxx:288
 AliTPCRF1D.cxx:289
 AliTPCRF1D.cxx:290
 AliTPCRF1D.cxx:291
 AliTPCRF1D.cxx:292
 AliTPCRF1D.cxx:293
 AliTPCRF1D.cxx:294
 AliTPCRF1D.cxx:295
 AliTPCRF1D.cxx:296
 AliTPCRF1D.cxx:297
 AliTPCRF1D.cxx:298
 AliTPCRF1D.cxx:299
 AliTPCRF1D.cxx:300
 AliTPCRF1D.cxx:301
 AliTPCRF1D.cxx:302
 AliTPCRF1D.cxx:303
 AliTPCRF1D.cxx:304
 AliTPCRF1D.cxx:305
 AliTPCRF1D.cxx:306
 AliTPCRF1D.cxx:307
 AliTPCRF1D.cxx:308
 AliTPCRF1D.cxx:309
 AliTPCRF1D.cxx:310
 AliTPCRF1D.cxx:311
 AliTPCRF1D.cxx:312
 AliTPCRF1D.cxx:313
 AliTPCRF1D.cxx:314
 AliTPCRF1D.cxx:315
 AliTPCRF1D.cxx:316
 AliTPCRF1D.cxx:317
 AliTPCRF1D.cxx:318
 AliTPCRF1D.cxx:319
 AliTPCRF1D.cxx:320
 AliTPCRF1D.cxx:321
 AliTPCRF1D.cxx:322
 AliTPCRF1D.cxx:323
 AliTPCRF1D.cxx:324
 AliTPCRF1D.cxx:325
 AliTPCRF1D.cxx:326
 AliTPCRF1D.cxx:327
 AliTPCRF1D.cxx:328
 AliTPCRF1D.cxx:329
 AliTPCRF1D.cxx:330
 AliTPCRF1D.cxx:331
 AliTPCRF1D.cxx:332
 AliTPCRF1D.cxx:333
 AliTPCRF1D.cxx:334
 AliTPCRF1D.cxx:335
 AliTPCRF1D.cxx:336
 AliTPCRF1D.cxx:337
 AliTPCRF1D.cxx:338
 AliTPCRF1D.cxx:339
 AliTPCRF1D.cxx:340
 AliTPCRF1D.cxx:341
 AliTPCRF1D.cxx:342
 AliTPCRF1D.cxx:343
 AliTPCRF1D.cxx:344
 AliTPCRF1D.cxx:345
 AliTPCRF1D.cxx:346
 AliTPCRF1D.cxx:347
 AliTPCRF1D.cxx:348
 AliTPCRF1D.cxx:349
 AliTPCRF1D.cxx:350
 AliTPCRF1D.cxx:351
 AliTPCRF1D.cxx:352
 AliTPCRF1D.cxx:353
 AliTPCRF1D.cxx:354
 AliTPCRF1D.cxx:355
 AliTPCRF1D.cxx:356
 AliTPCRF1D.cxx:357
 AliTPCRF1D.cxx:358
 AliTPCRF1D.cxx:359
 AliTPCRF1D.cxx:360
 AliTPCRF1D.cxx:361
 AliTPCRF1D.cxx:362
 AliTPCRF1D.cxx:363
 AliTPCRF1D.cxx:364
 AliTPCRF1D.cxx:365
 AliTPCRF1D.cxx:366
 AliTPCRF1D.cxx:367
 AliTPCRF1D.cxx:368
 AliTPCRF1D.cxx:369
 AliTPCRF1D.cxx:370
 AliTPCRF1D.cxx:371
 AliTPCRF1D.cxx:372
 AliTPCRF1D.cxx:373
 AliTPCRF1D.cxx:374
 AliTPCRF1D.cxx:375
 AliTPCRF1D.cxx:376
 AliTPCRF1D.cxx:377
 AliTPCRF1D.cxx:378
 AliTPCRF1D.cxx:379
 AliTPCRF1D.cxx:380
 AliTPCRF1D.cxx:381
 AliTPCRF1D.cxx:382
 AliTPCRF1D.cxx:383
 AliTPCRF1D.cxx:384
 AliTPCRF1D.cxx:385
 AliTPCRF1D.cxx:386
 AliTPCRF1D.cxx:387
 AliTPCRF1D.cxx:388
 AliTPCRF1D.cxx:389
 AliTPCRF1D.cxx:390
 AliTPCRF1D.cxx:391
 AliTPCRF1D.cxx:392
 AliTPCRF1D.cxx:393
 AliTPCRF1D.cxx:394
 AliTPCRF1D.cxx:395
 AliTPCRF1D.cxx:396
 AliTPCRF1D.cxx:397
 AliTPCRF1D.cxx:398
 AliTPCRF1D.cxx:399
 AliTPCRF1D.cxx:400
 AliTPCRF1D.cxx:401
 AliTPCRF1D.cxx:402
 AliTPCRF1D.cxx:403
 AliTPCRF1D.cxx:404