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.                  *
 **************************************************************************/

////
// This implementation AliTPCExB is using an "exact" calculation of the ExB
// effect. That is, it uses the drift ODE to calculate the distortion
// without any further assumption.
// Due to the numerical integration of the ODE, there are some numerical
// uncertencies involed.
////

#include "TMath.h"
#include "TTreeStream.h"
#include "AliMagF.h"
#include "AliTPCExBExact.h"

ClassImp(AliTPCExBExact)

const Double_t AliTPCExBExact::fgkEM=1.602176487e-19/9.10938215e-31;
const Double_t AliTPCExBExact::fgkDriftField=-40.e3;

AliTPCExBExact::AliTPCExBExact()
  : fDriftVelocity(0),
    //fkMap(0),
    fkField(0),fkN(0),
    fkNX(0),fkNY(0),fkNZ(0),
    fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
    fkZMin(-250.),fkZMax(250.),
    fkNLook(0),fkLook(0) {
  //
  // purely for I/O
  //
}

AliTPCExBExact::AliTPCExBExact(const AliMagF *bField,
			       Double_t driftVelocity,
			       Int_t nx,Int_t ny,Int_t nz,Int_t n)
  : fDriftVelocity(driftVelocity),
    //fkMap(0),
    fkField(bField),fkN(n),
    fkNX(nx),fkNY(ny),fkNZ(nz),
    fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
    fkZMin(-250.),fkZMax(250.),
    fkNLook(0),fkLook(0) {
  //
  // The constructor. One has to supply a magnetic field and an (initial)
  // drift velocity. Since some kind of lookuptable is created the
  // number of its meshpoints can be supplied.
  // n sets the number of integration steps to be used when integrating
  // over the full drift length.
  //
  CreateLookupTable();
}

/*
AliTPCExBExact::AliTPCExBExact(const AliFieldMap *bFieldMap,
			       Double_t driftVelocity,Int_t n) 
  : fDriftVelocity(driftVelocity),
    fkMap(bFieldMap),fkField(0),fkN(n), 
    fkNX(0),fkNY(0),fkNZ(0),
    fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
    fkZMin(-250.),fkZMax(250.),
    fkNLook(0),fkLook(0) {
  //
  // The constructor. One has to supply a field map and an (initial)
  // drift velocity.
  // n sets the number of integration steps to be used when integrating
  // over the full drift length.
  //

  fkXMin=bFieldMap->Xmin()
    -TMath::Ceil( (bFieldMap->Xmin()+250.0)/bFieldMap->DelX())
    *bFieldMap->DelX();
  fkXMax=bFieldMap->Xmax()
    -TMath::Floor((bFieldMap->Xmax()-250.0)/bFieldMap->DelX())
    *bFieldMap->DelX();
  fkYMin=bFieldMap->Ymin()
    -TMath::Ceil( (bFieldMap->Ymin()+250.0)/bFieldMap->DelY())
    *bFieldMap->DelY();
  fkYMax=bFieldMap->Ymax()
    -TMath::Floor((bFieldMap->Ymax()-250.0)/bFieldMap->DelY())
    *bFieldMap->DelY();
  fkZMax=bFieldMap->Zmax()
    -TMath::Floor((bFieldMap->Zmax()-250.0)/bFieldMap->DelZ())
    *bFieldMap->DelZ();
  fkZMax=TMath::Max(0.,fkZMax); // I really hope that this is unnecessary!

  fkNX=static_cast<Int_t>((fkXMax-fkXMin)/bFieldMap->DelX()+1.1);
  fkNY=static_cast<Int_t>((fkYMax-fkYMin)/bFieldMap->DelY()+1.1);
  fkNZ=static_cast<Int_t>((fkZMax-fkZMin)/bFieldMap->DelZ()+1.1);

  CreateLookupTable();
}
*/

AliTPCExBExact::~AliTPCExBExact() {
  //
  // destruct the poor object.
  //
  delete[] fkLook;
}

void AliTPCExBExact::Correct(const Double_t *position, Double_t *corrected) {
  //
  // correct for the distortion
  //
  Double_t x=(position[0]-fkXMin)/(fkXMax-fkXMin)*(fkNX-1);
  Int_t xi1=static_cast<Int_t>(x);
  xi1=TMath::Max(TMath::Min(xi1,fkNX-2),0);
  Int_t xi2=xi1+1;
  Double_t dx=(x-xi1);
  Double_t dx1=(xi2-x);
  
  Double_t y=(position[1]-fkYMin)/(fkYMax-fkYMin)*(fkNY-1);
  Int_t yi1=static_cast<Int_t>(y);
  yi1=TMath::Max(TMath::Min(yi1,fkNY-2),0);
  Int_t yi2=yi1+1;
  Double_t dy=(y-yi1);
  Double_t dy1=(yi2-y);
  
  Double_t z=position[2]/fkZMax*(fkNZ-1);
  Int_t side;
  if (z>0) {
    side=1;
  }
  else {
    z=-z;
    side=0;
  }
  Int_t zi1=static_cast<Int_t>(z);
  zi1=TMath::Max(TMath::Min(zi1,fkNZ-2),0);
  Int_t zi2=zi1+1;
  Double_t dz=(z-zi1);
  Double_t dz1=(zi2-z);
  
  for (int i=0;i<3;++i)
    corrected[i]
      =fkLook[(((xi1*fkNY+yi1)*fkNZ+zi1)*2+side)*3+i]*dx1*dy1*dz1
      +fkLook[(((xi1*fkNY+yi1)*fkNZ+zi2)*2+side)*3+i]*dx1*dy1*dz
      +fkLook[(((xi1*fkNY+yi2)*fkNZ+zi1)*2+side)*3+i]*dx1*dy *dz1
      +fkLook[(((xi1*fkNY+yi2)*fkNZ+zi2)*2+side)*3+i]*dx1*dy *dz
      +fkLook[(((xi2*fkNY+yi2)*fkNZ+zi1)*2+side)*3+i]*dx *dy *dz1
      +fkLook[(((xi2*fkNY+yi2)*fkNZ+zi2)*2+side)*3+i]*dx *dy *dz
      +fkLook[(((xi2*fkNY+yi1)*fkNZ+zi1)*2+side)*3+i]*dx *dy1*dz1
      +fkLook[(((xi2*fkNY+yi1)*fkNZ+zi2)*2+side)*3+i]*dx *dy1*dz ;
  //    corrected[2]=position[2];
}

/*
void AliTPCExBExact::TestThisBeautifulObject(const AliFieldMap *bFieldMap,
					     const char* fileName) {
  //
  // Have a look at the common part "TestThisBeautifulObjectGeneric".
  //
  fkMap=bFieldMap;
  fkField=0;
  TestThisBeautifulObjectGeneric(fileName);
}
*/

void AliTPCExBExact::TestThisBeautifulObject(const AliMagF *bField,
					     const char* fileName) {
  //
  // Have a look at the common part "TestThisBeautifulObjectGeneric".
  //
  fkField=bField;
  //fkMap=0;
  TestThisBeautifulObjectGeneric(fileName);
}

void AliTPCExBExact::TestThisBeautifulObjectGeneric(const char* fileName) {
  //
  // Well, as the name sais... it tests the object.
  //
  TTreeSRedirector ts(fileName);
  Double_t x[3];
  for (x[0]=-250.;x[0]<=250.;x[0]+=10.)
    for (x[1]=-250.;x[1]<=250.;x[1]+=10.)
      for (x[2]=-250.;x[2]<=250.;x[2]+=10.) {
	Double_t d[3];
	Double_t dnl[3];
	Correct(x,d);
	CalculateDistortion(x,dnl);
	Double_t r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
	Double_t rd=TMath::Sqrt(d[0]*d[0]+d[1]*d[1]);
	Double_t dr=r-rd;
	Double_t phi=TMath::ATan2(x[0],x[1]);
	Double_t phid=TMath::ATan2(d[0],d[1]);
	Double_t dphi=phi-phid;
	if (dphi<0.) dphi+=TMath::TwoPi();
	if (dphi>TMath::Pi()) dphi=TMath::TwoPi()-dphi;
	Double_t drphi=r*dphi;
	Double_t dx=x[0]-d[0];
	Double_t dy=x[1]-d[1];
	Double_t dz=x[2]-d[2];
	Double_t dnlx=x[0]-dnl[0];
	Double_t dnly=x[1]-dnl[1];
	Double_t dnlz=x[2]-dnl[2];
	Double_t b[3];
	GetB(b,x);
	ts<<"positions"
	  <<"bx="<<b[0]
	  <<"by="<<b[1]
	  <<"bz="<<b[2]
	  <<"x0="<<x[0]
	  <<"x1="<<x[1]
	  <<"x2="<<x[2]
	  <<"dx="<<dx
	  <<"dy="<<dy
	  <<"dz="<<dz
	  <<"dnlx="<<dnlx
	  <<"dnly="<<dnly
	  <<"dnlz="<<dnlz
	  <<"r="<<r
	  <<"phi="<<phi
	  <<"dr="<<dr
	  <<"drphi="<<drphi
	  <<"\n";
      }
}

void AliTPCExBExact::CreateLookupTable() {
  //
  // Helper function to fill the lookup table.
  //
  fkNLook=fkNX*fkNY*fkNZ*2*3;
  fkLook=new Double_t[fkNLook];
  Double_t x[3];
  for (int i=0;i<fkNX;++i) {
    x[0]=fkXMin+(fkXMax-fkXMin)/(fkNX-1)*i;
    for (int j=0;j<fkNY;++j) {
      x[1]=fkYMin+(fkYMax-fkYMin)/(fkNY-1)*j;
      for (int k=0;k<fkNZ;++k) {
	x[2]=1.*fkZMax/(fkNZ-1)*k;
	x[2]=TMath::Max((Double_t)0.0001,x[2]); //ugly
	CalculateDistortion(x,&fkLook[(((i*fkNY+j)*fkNZ+k)*2+1)*3]);
	x[2]=-x[2];
	CalculateDistortion(x,&fkLook[(((i*fkNY+j)*fkNZ+k)*2+0)*3]);
      }
    }
  }
}

void AliTPCExBExact::GetE(Double_t *e,const Double_t *x) const {
  //
  // Helper function returning the E field in SI units (V/m).
  //
  e[0]=0.;
  e[1]=0.;
  e[2]=(x[2]<0.?-1.:1.)*fgkDriftField; // in V/m
}

void AliTPCExBExact::GetB(Double_t *b,const Double_t *x) const {
  //
  // Helper function returning the B field in SI units (T).
  //
  Double_t xm[3];
  // the beautiful m to cm (and the ugly "const_cast") and Double_t 
  // to Float_t read the NRs introduction!:
  for (int i=0;i<3;++i) xm[i]=x[i]*100.;
  Double_t bf[3];
  //if (fkMap!=0)
  //  fkMap->Field(xm,bf);
  //else
  ((AliMagF*)fkField)->Field(xm,bf);
  for (int i=0;i<3;++i) b[i]=bf[i]/10.;
}

void AliTPCExBExact::Motion(const Double_t *x,Double_t,
			    Double_t *dxdt) const {
  //
  // The differential equation of motion of the electrons.
  //
  Double_t tau=fDriftVelocity/fgkDriftField/fgkEM;
  Double_t tau2=tau*tau;
  Double_t e[3];
  Double_t b[3];
  GetE(e,x);
  GetB(b,x);
  Double_t wx=fgkEM*b[0];
  Double_t wy=fgkEM*b[1];
  Double_t wz=fgkEM*b[2];
  Double_t ex=fgkEM*e[0];
  Double_t ey=fgkEM*e[1];
  Double_t ez=fgkEM*e[2];
  Double_t w2=(wx*wx+wy*wy+wz*wz);
  dxdt[0]=(1.+wx*wx*tau2)*ex+(wz*tau+wx*wy*tau2)*ey+(-wy*tau+wx*wz*tau2)*ez;
  dxdt[1]=(-wz*tau+wx*wy*tau2)*ex+(1.+wy*wy*tau2)*ey+(wx*tau+wy*wz*tau2)*ez;
  dxdt[2]=(wy*tau+wx*wz*tau2)*ex+(-wx*tau+wy*wz*tau2)*ey+(1.+wz*wz*tau2)*ez;
  Double_t fac=tau/(1.+w2*tau2);
  dxdt[0]*=fac;
  dxdt[1]*=fac;
  dxdt[2]*=fac;
}

void AliTPCExBExact::CalculateDistortion(const Double_t *x0,
					 Double_t *dist) const {
  //
  // Helper function that calculates one distortion by integration
  // (only used to fill the lookup table).
  //
  Double_t h=0.01*250./fDriftVelocity/fkN;
  Double_t t=0.;
  Double_t xt[3];
  Double_t xo[3];
  for (int i=0;i<3;++i)
    xo[i]=xt[i]=x0[i]*0.01;
  while (TMath::Abs(xt[2])<250.*0.01) {
    for (int i=0;i<3;++i)
      xo[i]=xt[i];
    DGLStep(xt,t,h);
    t+=h;
  }
  if (t!=0.) {
    Double_t p=((xt[2]<0.?-1.:1.)*250.*0.01-xo[2])/(xt[2]-xo[2]);
    dist[0]=(xo[0]+p*(xt[0]-xo[0]))*100.;
    dist[1]=(xo[1]+p*(xt[1]-xo[1]))*100.;
    //    dist[2]=(xo[2]+p*(xt[2]-xo[2]))*100.;
    dist[2]=(x0[2]>0.?-1:1.)*(t-h+p*h)*fDriftVelocity*100.;
    dist[2]+=(x0[2]<0.?-1:1.)*250.;
  }
  else {
    dist[0]=x0[0];
    dist[1]=x0[1];
    dist[2]=x0[2];
  }
  // reverse the distortion, i.e. get the correction
  dist[0]=x0[0]-(dist[0]-x0[0]);
  dist[1]=x0[1]-(dist[1]-x0[1]);
}

void AliTPCExBExact::DGLStep(Double_t *x,Double_t t,Double_t h) const {
  //
  // An elementary integration step.
  // (simple Euler Method)
  //
  Double_t dxdt[3];
  Motion(x,t,dxdt);
  for (int i=0;i<3;++i)
    x[i]+=h*dxdt[i];

  /* suggestions about how to write it this way are welcome!
     void DGLStep(void (*f)(const Double_t *x,Double_t t,Double_t *dxdt),
                   Double_t *x,Double_t t,Double_t h,Int_t n) const;
  */

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