ROOT logo
///////////////////////////////////////////////////////////////////////////
//                                                                       //
// AliFemtoCoulomb: This is a Coulomb correction class which             //
//  1. Reads in the dat from a file                                      //  
//  2. Performs a linear interpolation in R and creates any array of     //
//     interpolations                                                    //
//  3. Interpolates in eta and returns the Coulomb correction to user    //
//                                                                       //
///////////////////////////////////////////////////////////////////////////

#include "AliFemtoCoulomb.h"
//#include "Stiostream.h"
#include <stdio.h>
#include <cassert>
//#include "PhysicalConstants.h"
#define fine_structure_const 0.00729735

#ifdef __ROOT__
ClassImp(AliFemtoCoulomb)
#endif

AliFemtoCoulomb::AliFemtoCoulomb() :
  fFile(""),
  fRadius(-1.0),
  fZ1Z2(1.0),
  fNLines(0)
{
  // Default constructor
  fFile = "/afs/rhic/star/hbt/coul/AliFemtoCorrectionFiles/correctionpp.dat";
  if (!fFile) {
    cout << " No file, dummy!" << endl;
    assert(0);
  }
  cout << "You have 1 default Coulomb correction!" << endl;

  for (int ie=0; ie<1000; ie++) {
    fEta[ie] = 0.0;
    fCoulomb[ie] = 0.0;
  }
}

AliFemtoCoulomb::AliFemtoCoulomb(const AliFemtoCoulomb& aCoul) :
  fFile(aCoul.fFile),
  fRadius(aCoul.fRadius),
  fZ1Z2(aCoul.fZ1Z2),
  fNLines(0)
{
  // copy constructor
  CreateLookupTable(fRadius);
}

AliFemtoCoulomb::AliFemtoCoulomb(const char* readFile, const double& radius, const double& charge) :
  fFile(readFile),
  fRadius(radius),
  fZ1Z2(0),
  fNLines(0)
{
  // constructor with explicit filename
  fFile = readFile;
  fRadius = radius;
  CreateLookupTable(fRadius);
  fZ1Z2 = charge;
  cout << "You have 1 Coulomb correction!" << endl;
}

AliFemtoCoulomb::~AliFemtoCoulomb() {
  // destructor
}

AliFemtoCoulomb& AliFemtoCoulomb::operator=(const AliFemtoCoulomb& aCoul)
{
  // assignment operator
  if (this == &aCoul)
    return *this;

  fFile = aCoul.fFile;
  fRadius = aCoul.fRadius;
  fZ1Z2 = aCoul.fZ1Z2;

  CreateLookupTable(fRadius);
  
  return *this;
}


void AliFemtoCoulomb::SetRadius(const double& radius) {
  // set the coulomb radius
  cout << " AliFemtoCoulomb::setRadius() " << endl;
  fRadius = radius;
  CreateLookupTable(fRadius);
}

double AliFemtoCoulomb::GetRadius() const {
  // return coulomb radius
  return (fRadius);
}

void AliFemtoCoulomb::SetFile(const char* readFile) {
  // set the filename with coulomb calculations
  cout << " AliFemtoCoulomb::SetFile() " << endl;
  fFile = readFile;
  // Create new lookup table since file has changed
  if (fRadius>0.0) {
    CreateLookupTable(fRadius);
  }
}

void AliFemtoCoulomb::SetChargeProduct(const double& charge) {
  // set pair charge
  cout << " AliFemtoCoulomb::SetChargeProduct() " << endl;
  if ( fZ1Z2!=charge ) { 
    fZ1Z2 = charge;
    if ( fZ1Z2>0 ) {
      fFile = "/afs/rhic/star/hbt/coul/AliFemtoCorrectionFiles/correctionpp.dat";
    }
    else {
      fFile = "/afs/rhic/star/hbt/coul/AliFemtoCorrectionFiles/correctionpm.dat";
    }
    CreateLookupTable(fRadius);
  }
}

void AliFemtoCoulomb::CreateLookupTable(const double& radius) {
  // Read radii from fFile
  // Create array(pair) of linear interpolation between radii
  cout << " AliFemtoCoulomb::CreateLookupTable() " << endl;

  if (radius<0.0) {
    cout << " AliFemtoCoulomb::CreateLookupTable -> NEGATIVE RADIUS " << endl;
    cout << "  call AliFemtoCoulomb::SetRadius(r) with positive r " << endl;
    cerr << " AliFemtoCoulomb::CreateLookupTable -> NEGATIVE RADIUS " << endl;
    cerr << "  call AliFemtoCoulomb::SetRadius(r) with positive r " << endl;
    assert(0);
  }
  ifstream mystream(fFile);
  if (!mystream) {
    cout << "Could not open file" << endl;
    assert(0);
  }
  else {
    cout << "Input correction file opened" << endl;
  }

  static char tempstring[2001];
  static float radii[2000];
  static int tNRadii = 0;
  tNRadii = 0;
  if (!mystream.getline(tempstring,2000)) {
    cout << "Could not read radii from file" << endl;
    assert(0);
  }
  for (unsigned int ii=0; ii<strlen(tempstring); ii++) {
    while (tempstring[ii]==' ') ii++;
    sscanf(&tempstring[ii++],"%f",&radii[++tNRadii]);
    while ( tempstring[ii]!=' ' && (ii)<strlen(tempstring) )ii++;
  }
  cout << " Read " << tNRadii << " radii from file" << endl;

  static double tLowRadius = -1.0;
  static double tHighRadius = -1.0;
  static int tLowIndex = 0;
  tLowRadius = -1.0;
  tHighRadius = -1.0;
  tLowIndex = 0;
  for(int iii=1; iii<=tNRadii-1; iii++) { // Loop to one less than #radii
    if ( radius >= radii[iii] && radius <= radii[iii+1] ) {
      tLowRadius = radii[iii];
      tHighRadius = radii[iii+1];
      tLowIndex = iii;
    }
  }
  if ( (tLowRadius < 0.0) || (tHighRadius < 0.0) ) {
    cout << "AliFemtoCoulomb::CreateLookupTable --> Problem interpolating radius" << endl;
    cout << "  Check range of radii in lookup file...." << endl;
    cerr << "AliFemtoCoulomb::CreateLookupTable --> Problem interpolating radius" << endl;
    cerr << "  Check range of radii in lookup file...." << endl;
    assert(0);
  }

  static double corr[100];           // array of corrections ... must be > tNRadii
  fNLines = 0;
  static double tempEta = 0;
  tempEta = 0;
  while (mystream >> tempEta) {
    for (int i=1; i<=tNRadii; i++) {
      mystream >> corr[i];
    }
    static double tLowCoulomb = 0;
    static double tHighCoulomb = 0;
    static double nCorr = 0;
    tLowCoulomb = corr[tLowIndex];
    tHighCoulomb = corr[tLowIndex+1];
    nCorr = ( (radius-tLowRadius)*tHighCoulomb+(tHighRadius-radius)*tLowCoulomb )/(tHighRadius-tLowRadius);
      fEta[fNLines] = tempEta;     // Eta
      fCoulomb[fNLines] = nCorr;   // Interpolated Coulomb correction for radius
      fNLines++;
  }
  mystream.close();
  cout << "Lookup Table is created with " << fNLines << " points" << endl;
}

double AliFemtoCoulomb::CoulombCorrect(const double& eta) {
  // Interpolates in eta
  if (fRadius < 0.0) {
    cout << "AliFemtoCoulomb::CoulombCorrect(eta) --> Trying to correct for negative radius!" << endl;
    cerr << "AliFemtoCoulomb::CoulombCorrect(eta) --> Trying to correct for negative radius!" << endl;
    assert(0);
  }
  static int middle=0;
  middle=int( (fNLines-1)/2 );
  if (eta*fEta[middle]<0.0) {
    cout << "AliFemtoCoulomb::CoulombCorrect(eta) --> eta: " << eta << " has wrong sign for data file! " << endl;
    cerr << "AliFemtoCoulomb::CoulombCorrect(eta) --> eta: " << eta << " has wrong sign for data file! " << endl;
    assert(0);
  }

  static double tCorr = 0;
  tCorr = -1.0;
  
  if ( (eta>fEta[0]) && (fEta[0]>0.0) ) {
    tCorr = fCoulomb[0];
    return (tCorr);
  }
  if ( (eta<fEta[fNLines-1]) && (fEta[fNLines-1]<0.0) ) {
    tCorr = fCoulomb[fNLines-1];
    return (tCorr);
  }
  // This is a binary search for the bracketing pair of data points
  static int high = 0;
  static int low = 0;
  static int width = 0;
  high = fNLines-1;
  low = 0;
  width = high-low;
  middle = int(width/2.0); // Was instantiated above
  while (middle > 0) {
    if (fEta[low+middle] < eta) {
      // eta is in the 1st half
      high-=middle;
      width = high-low;
      middle = int(width/2.0);
    }
    else {
      // eta is in the 2nd half
      low+=middle;
      width = high-low;
      middle = int(width/2.0);
    }
  }
  // Make sure we found the right one
  if ( (fEta[low] >= eta) && (eta >= fEta[low+1]) ) {
    static double tLowEta = 0;
    static double tHighEta = 0;    
    static double tLowCoulomb = 0;
    static double tHighCoulomb = 0;
    tLowEta = fEta[low];
    tHighEta = fEta[low+1];    
    tLowCoulomb = fCoulomb[low];
    tHighCoulomb = fCoulomb[low+1];
    //      cout << tLowEta << " *** Eta *** " << tHighEta << endl;
    //      cout << tLowCoulomb << " *** Coulomb *** " << tHighCoulomb << endl;
    tCorr = ( (eta-tLowEta)*tHighCoulomb+(tHighEta-eta)*tLowCoulomb )/(tHighEta-tLowEta);
  }
  if (tCorr<0.0) {
    cout << "AliFemtoCoulomb::CoulombCorrect(eta) --> No correction" << endl;
    cout << "  Check range of eta in file: Input eta  " << eta << endl;
    cerr << "AliFemtoCoulomb::CoulombCorrect(eta) --> No correction" << endl;
    cerr << "  Check range of eta in file: Input eta  " << eta << endl;
    assert(0);
  } 
  return (tCorr);

}

double AliFemtoCoulomb::CoulombCorrect(const double& eta,
				const double& radius) {
  // Checks radii ... input radius and fRadius
  // Calls createLookupTable if neccessary
  // Interpolate(linear) between etas in the created lookup table

  if (radius < 0.0) {
    if (fRadius < 0.0) {
      // Both radii are negative
      cout << "AliFemtoCoulomb::CoulombCorrect(eta,r) --> input and member radii are negative!" << endl;
      cerr << "AliFemtoCoulomb::CoulombCorrect(eta,r) --> input and member radii are negative!" << endl;
      assert(0);
    }
  }
  else {
    // radius > 0.0
    if (radius == fRadius) {
      // Both radii are positive and equal
      //      cout << "Radii are the same!!!" << endl;
    }
    else {
      // Both radii are positive but not equal
      fRadius = radius;
      CreateLookupTable(fRadius);
    }
  }

  // Interpolate in eta
  return ( CoulombCorrect(eta) );
}

double AliFemtoCoulomb::CoulombCorrect(const AliFemtoPair* pair) {
  return ( CoulombCorrect( Eta(pair) ) );;
}

double AliFemtoCoulomb::CoulombCorrect(const AliFemtoPair* pair, const double& radius) {
  return ( CoulombCorrect( Eta(pair),radius ) );
}

double AliFemtoCoulomb::Eta(const AliFemtoPair* pair) {
  // calculate eta
  static double px1,py1,pz1,px2,py2,pz2;
  static double px1new,py1new,pz1new;
  static double px2new,py2new,pz2new;
  static double vx1cms,vy1cms,vz1cms;
  static double vx2cms,vy2cms,vz2cms;
  static double tVcmsX,tVcmsY,tVcmsZ;
  static double dv = 0.0;
  static double e1,e2,e1new,e2new;
  static double psi,theta;
  static double beta,gamma;
  static double tVcmsXnew;

  px1 = pair->Track1()->FourMomentum().px();
  py1 = pair->Track1()->FourMomentum().py();
  pz1 = pair->Track1()->FourMomentum().pz();
  e1 = pair->Track1()->FourMomentum().e();
  px2 = pair->Track2()->FourMomentum().px();
  py2 = pair->Track2()->FourMomentum().py();
  pz2 = pair->Track2()->FourMomentum().pz();
  e2 = pair->Track2()->FourMomentum().e();
  
  tVcmsX = ( px1+px2 )/( e1+e2 );
  tVcmsY = ( py1+py2 )/( e1+e2 );
  tVcmsZ = ( pz1+pz2 )/( e1+e2 );
  // Rotate tVcms to x-direction
  psi = atan(tVcmsY/tVcmsX);
  tVcmsXnew = tVcmsX*cos(psi)+tVcmsY*sin(psi);
  tVcmsX = tVcmsXnew;
  theta = atan(tVcmsZ/tVcmsX);
  tVcmsXnew = tVcmsX*cos(theta)+tVcmsZ*sin(theta);
  tVcmsX = tVcmsXnew;
  // Gamma and Beta
  beta = tVcmsX;
  gamma = 1.0/::sqrt( 1.0-beta*beta );

  // Rotate p1 and p2 to new frame
  px1new = px1*cos(psi)+py1*sin(psi);
  py1new = -px1*sin(psi)+py1*cos(psi);
  px1 = px1new;
  px1new = px1*cos(theta)+pz1*sin(theta);
  pz1new = -px1*sin(theta)+pz1*cos(theta);
  px1 = px1new;
  py1 = py1new;
  pz1 = pz1new;

  px2new = px2*cos(psi)+py2*sin(psi);
  py2new = -px2*sin(psi)+py2*cos(psi);
  px2 = px2new;
  px2new = px2*cos(theta)+pz2*sin(theta);
  pz2new = -px2*sin(theta)+pz2*cos(theta);
  px2 = px2new;
  py2 = py2new;
  pz2 = pz2new;

  // Lorentz transform the x component and energy
  e1new = gamma*e1 - gamma*beta*px1;
  px1new = -gamma*beta*e1 + gamma*px1;
  e2new = gamma*e2 - gamma*beta*px2;
  px2new = -gamma*beta*e2 + gamma*px2;
  px1 = px1new;
  px2 = px2new;

  // New velocities
  vx1cms = px1/e1new;
  vy1cms = py1/e1new;
  vz1cms = pz1/e1new;
  vx2cms = px2/e2new;
  vy2cms = py2/e2new;
  vz2cms = pz2/e2new;

  // Velocity difference in CMS frame
  dv = ::sqrt( (vx1cms-vx2cms)*(vx1cms-vx2cms) +
	     (vy1cms-vy2cms)*(vy1cms-vy2cms) +
	     (vz1cms-vz2cms)*(vz1cms-vz2cms) );

  return ( fZ1Z2*fine_structure_const/(dv) );
}

TH1D* AliFemtoCoulomb::CorrectionHistogram(const double& mass1, const double& mass2, const int& nBins, 
						const double& low, const double& high) {
  // return correction histogram

  if ( mass1!=mass2 ) {
    cout << "Masses not equal ... try again.  No histogram created." << endl;
    assert(0);
  }
  TH1D* correction = new TH1D("correction","Coulomb correction",nBins,low,high);
  const double kReducedMass = mass1*mass2/(mass1+mass2);
  double qInv = low;
  //double dQinv = (high-low)/( (double)nBins );
  double eta;
  for (int ii=1; ii<=nBins; ii++) 
    {
      qInv = correction->GetBinCenter(ii);
      eta = 2.0*fZ1Z2*kReducedMass*fine_structure_const/( qInv );
      CoulombCorrect( eta );
      correction->Fill( qInv, CoulombCorrect(eta,fRadius) );
    }

  return (correction);
}

#ifdef __ROOT__
TH1D* AliFemtoCoulomb::CorrectionHistogram(const TH1D* histo, const double mass) {
  // return correction histogram - 1D case
  TH1D* correction = (TH1D*) ((TH1D*)histo)->Clone();
  correction->Reset();
  correction->SetDirectory(0);
  int    nBins = correction->GetXaxis()->GetNbins();
  const double kReducedMass = 0.5*mass;
  double qInv;
  double eta;
  for (int ii=1; ii<=nBins; ii++) 
    {
      qInv = correction->GetBinCenter(ii);
      eta = 2.0*fZ1Z2*kReducedMass*fine_structure_const/( qInv );
      correction->Fill( qInv, CoulombCorrect(eta,fRadius) );
    }

  return (correction);
}

TH3D* AliFemtoCoulomb::CorrectionHistogram(const TH3D* histo, const double mass) {
  // return correction histogram - 3D case
  TH3D* correction = (TH3D*) ((TH3D*)histo)->Clone();
  correction->Reset();
  correction->SetDirectory(0);
  int    nBinsX = correction->GetXaxis()->GetNbins();
  int    nBinsY = correction->GetYaxis()->GetNbins();
  int    nBinsZ = correction->GetZaxis()->GetNbins();
  const double kReducedMass = 0.5*mass;
  double eta;
  double qInv;
  int binNumber;
  for (int ii=1; ii<=nBinsX; ii++) { 
    for (int iii=1; iii<=nBinsY; iii++) {
      for (int iv=1; iv<=nBinsZ; iv++) {
	binNumber = histo->GetBin(ii,iii,iv);
	qInv = histo->GetBinContent(binNumber);
	eta = 2.0*fZ1Z2*kReducedMass*fine_structure_const/( qInv );
	correction->SetBinContent(binNumber, CoulombCorrect(eta,fRadius) );
      }
    }
  }
  return (correction);
}
#endif

double AliFemtoCoulomb::CoulombCorrect(const double& mass, const double& charge,
				    const double& radius, const double& qInv) {
  // return correction factor
  fRadius = radius;
  fZ1Z2 = charge;
  const double kReducedMass = 0.5*mass; // must be same mass particles
  double eta = 2.0*fZ1Z2*kReducedMass*fine_structure_const/( qInv );
  return ( CoulombCorrect(eta,fRadius) );
}
 AliFemtoCoulomb.cxx:1
 AliFemtoCoulomb.cxx:2
 AliFemtoCoulomb.cxx:3
 AliFemtoCoulomb.cxx:4
 AliFemtoCoulomb.cxx:5
 AliFemtoCoulomb.cxx:6
 AliFemtoCoulomb.cxx:7
 AliFemtoCoulomb.cxx:8
 AliFemtoCoulomb.cxx:9
 AliFemtoCoulomb.cxx:10
 AliFemtoCoulomb.cxx:11
 AliFemtoCoulomb.cxx:12
 AliFemtoCoulomb.cxx:13
 AliFemtoCoulomb.cxx:14
 AliFemtoCoulomb.cxx:15
 AliFemtoCoulomb.cxx:16
 AliFemtoCoulomb.cxx:17
 AliFemtoCoulomb.cxx:18
 AliFemtoCoulomb.cxx:19
 AliFemtoCoulomb.cxx:20
 AliFemtoCoulomb.cxx:21
 AliFemtoCoulomb.cxx:22
 AliFemtoCoulomb.cxx:23
 AliFemtoCoulomb.cxx:24
 AliFemtoCoulomb.cxx:25
 AliFemtoCoulomb.cxx:26
 AliFemtoCoulomb.cxx:27
 AliFemtoCoulomb.cxx:28
 AliFemtoCoulomb.cxx:29
 AliFemtoCoulomb.cxx:30
 AliFemtoCoulomb.cxx:31
 AliFemtoCoulomb.cxx:32
 AliFemtoCoulomb.cxx:33
 AliFemtoCoulomb.cxx:34
 AliFemtoCoulomb.cxx:35
 AliFemtoCoulomb.cxx:36
 AliFemtoCoulomb.cxx:37
 AliFemtoCoulomb.cxx:38
 AliFemtoCoulomb.cxx:39
 AliFemtoCoulomb.cxx:40
 AliFemtoCoulomb.cxx:41
 AliFemtoCoulomb.cxx:42
 AliFemtoCoulomb.cxx:43
 AliFemtoCoulomb.cxx:44
 AliFemtoCoulomb.cxx:45
 AliFemtoCoulomb.cxx:46
 AliFemtoCoulomb.cxx:47
 AliFemtoCoulomb.cxx:48
 AliFemtoCoulomb.cxx:49
 AliFemtoCoulomb.cxx:50
 AliFemtoCoulomb.cxx:51
 AliFemtoCoulomb.cxx:52
 AliFemtoCoulomb.cxx:53
 AliFemtoCoulomb.cxx:54
 AliFemtoCoulomb.cxx:55
 AliFemtoCoulomb.cxx:56
 AliFemtoCoulomb.cxx:57
 AliFemtoCoulomb.cxx:58
 AliFemtoCoulomb.cxx:59
 AliFemtoCoulomb.cxx:60
 AliFemtoCoulomb.cxx:61
 AliFemtoCoulomb.cxx:62
 AliFemtoCoulomb.cxx:63
 AliFemtoCoulomb.cxx:64
 AliFemtoCoulomb.cxx:65
 AliFemtoCoulomb.cxx:66
 AliFemtoCoulomb.cxx:67
 AliFemtoCoulomb.cxx:68
 AliFemtoCoulomb.cxx:69
 AliFemtoCoulomb.cxx:70
 AliFemtoCoulomb.cxx:71
 AliFemtoCoulomb.cxx:72
 AliFemtoCoulomb.cxx:73
 AliFemtoCoulomb.cxx:74
 AliFemtoCoulomb.cxx:75
 AliFemtoCoulomb.cxx:76
 AliFemtoCoulomb.cxx:77
 AliFemtoCoulomb.cxx:78
 AliFemtoCoulomb.cxx:79
 AliFemtoCoulomb.cxx:80
 AliFemtoCoulomb.cxx:81
 AliFemtoCoulomb.cxx:82
 AliFemtoCoulomb.cxx:83
 AliFemtoCoulomb.cxx:84
 AliFemtoCoulomb.cxx:85
 AliFemtoCoulomb.cxx:86
 AliFemtoCoulomb.cxx:87
 AliFemtoCoulomb.cxx:88
 AliFemtoCoulomb.cxx:89
 AliFemtoCoulomb.cxx:90
 AliFemtoCoulomb.cxx:91
 AliFemtoCoulomb.cxx:92
 AliFemtoCoulomb.cxx:93
 AliFemtoCoulomb.cxx:94
 AliFemtoCoulomb.cxx:95
 AliFemtoCoulomb.cxx:96
 AliFemtoCoulomb.cxx:97
 AliFemtoCoulomb.cxx:98
 AliFemtoCoulomb.cxx:99
 AliFemtoCoulomb.cxx:100
 AliFemtoCoulomb.cxx:101
 AliFemtoCoulomb.cxx:102
 AliFemtoCoulomb.cxx:103
 AliFemtoCoulomb.cxx:104
 AliFemtoCoulomb.cxx:105
 AliFemtoCoulomb.cxx:106
 AliFemtoCoulomb.cxx:107
 AliFemtoCoulomb.cxx:108
 AliFemtoCoulomb.cxx:109
 AliFemtoCoulomb.cxx:110
 AliFemtoCoulomb.cxx:111
 AliFemtoCoulomb.cxx:112
 AliFemtoCoulomb.cxx:113
 AliFemtoCoulomb.cxx:114
 AliFemtoCoulomb.cxx:115
 AliFemtoCoulomb.cxx:116
 AliFemtoCoulomb.cxx:117
 AliFemtoCoulomb.cxx:118
 AliFemtoCoulomb.cxx:119
 AliFemtoCoulomb.cxx:120
 AliFemtoCoulomb.cxx:121
 AliFemtoCoulomb.cxx:122
 AliFemtoCoulomb.cxx:123
 AliFemtoCoulomb.cxx:124
 AliFemtoCoulomb.cxx:125
 AliFemtoCoulomb.cxx:126
 AliFemtoCoulomb.cxx:127
 AliFemtoCoulomb.cxx:128
 AliFemtoCoulomb.cxx:129
 AliFemtoCoulomb.cxx:130
 AliFemtoCoulomb.cxx:131
 AliFemtoCoulomb.cxx:132
 AliFemtoCoulomb.cxx:133
 AliFemtoCoulomb.cxx:134
 AliFemtoCoulomb.cxx:135
 AliFemtoCoulomb.cxx:136
 AliFemtoCoulomb.cxx:137
 AliFemtoCoulomb.cxx:138
 AliFemtoCoulomb.cxx:139
 AliFemtoCoulomb.cxx:140
 AliFemtoCoulomb.cxx:141
 AliFemtoCoulomb.cxx:142
 AliFemtoCoulomb.cxx:143
 AliFemtoCoulomb.cxx:144
 AliFemtoCoulomb.cxx:145
 AliFemtoCoulomb.cxx:146
 AliFemtoCoulomb.cxx:147
 AliFemtoCoulomb.cxx:148
 AliFemtoCoulomb.cxx:149
 AliFemtoCoulomb.cxx:150
 AliFemtoCoulomb.cxx:151
 AliFemtoCoulomb.cxx:152
 AliFemtoCoulomb.cxx:153
 AliFemtoCoulomb.cxx:154
 AliFemtoCoulomb.cxx:155
 AliFemtoCoulomb.cxx:156
 AliFemtoCoulomb.cxx:157
 AliFemtoCoulomb.cxx:158
 AliFemtoCoulomb.cxx:159
 AliFemtoCoulomb.cxx:160
 AliFemtoCoulomb.cxx:161
 AliFemtoCoulomb.cxx:162
 AliFemtoCoulomb.cxx:163
 AliFemtoCoulomb.cxx:164
 AliFemtoCoulomb.cxx:165
 AliFemtoCoulomb.cxx:166
 AliFemtoCoulomb.cxx:167
 AliFemtoCoulomb.cxx:168
 AliFemtoCoulomb.cxx:169
 AliFemtoCoulomb.cxx:170
 AliFemtoCoulomb.cxx:171
 AliFemtoCoulomb.cxx:172
 AliFemtoCoulomb.cxx:173
 AliFemtoCoulomb.cxx:174
 AliFemtoCoulomb.cxx:175
 AliFemtoCoulomb.cxx:176
 AliFemtoCoulomb.cxx:177
 AliFemtoCoulomb.cxx:178
 AliFemtoCoulomb.cxx:179
 AliFemtoCoulomb.cxx:180
 AliFemtoCoulomb.cxx:181
 AliFemtoCoulomb.cxx:182
 AliFemtoCoulomb.cxx:183
 AliFemtoCoulomb.cxx:184
 AliFemtoCoulomb.cxx:185
 AliFemtoCoulomb.cxx:186
 AliFemtoCoulomb.cxx:187
 AliFemtoCoulomb.cxx:188
 AliFemtoCoulomb.cxx:189
 AliFemtoCoulomb.cxx:190
 AliFemtoCoulomb.cxx:191
 AliFemtoCoulomb.cxx:192
 AliFemtoCoulomb.cxx:193
 AliFemtoCoulomb.cxx:194
 AliFemtoCoulomb.cxx:195
 AliFemtoCoulomb.cxx:196
 AliFemtoCoulomb.cxx:197
 AliFemtoCoulomb.cxx:198
 AliFemtoCoulomb.cxx:199
 AliFemtoCoulomb.cxx:200
 AliFemtoCoulomb.cxx:201
 AliFemtoCoulomb.cxx:202
 AliFemtoCoulomb.cxx:203
 AliFemtoCoulomb.cxx:204
 AliFemtoCoulomb.cxx:205
 AliFemtoCoulomb.cxx:206
 AliFemtoCoulomb.cxx:207
 AliFemtoCoulomb.cxx:208
 AliFemtoCoulomb.cxx:209
 AliFemtoCoulomb.cxx:210
 AliFemtoCoulomb.cxx:211
 AliFemtoCoulomb.cxx:212
 AliFemtoCoulomb.cxx:213
 AliFemtoCoulomb.cxx:214
 AliFemtoCoulomb.cxx:215
 AliFemtoCoulomb.cxx:216
 AliFemtoCoulomb.cxx:217
 AliFemtoCoulomb.cxx:218
 AliFemtoCoulomb.cxx:219
 AliFemtoCoulomb.cxx:220
 AliFemtoCoulomb.cxx:221
 AliFemtoCoulomb.cxx:222
 AliFemtoCoulomb.cxx:223
 AliFemtoCoulomb.cxx:224
 AliFemtoCoulomb.cxx:225
 AliFemtoCoulomb.cxx:226
 AliFemtoCoulomb.cxx:227
 AliFemtoCoulomb.cxx:228
 AliFemtoCoulomb.cxx:229
 AliFemtoCoulomb.cxx:230
 AliFemtoCoulomb.cxx:231
 AliFemtoCoulomb.cxx:232
 AliFemtoCoulomb.cxx:233
 AliFemtoCoulomb.cxx:234
 AliFemtoCoulomb.cxx:235
 AliFemtoCoulomb.cxx:236
 AliFemtoCoulomb.cxx:237
 AliFemtoCoulomb.cxx:238
 AliFemtoCoulomb.cxx:239
 AliFemtoCoulomb.cxx:240
 AliFemtoCoulomb.cxx:241
 AliFemtoCoulomb.cxx:242
 AliFemtoCoulomb.cxx:243
 AliFemtoCoulomb.cxx:244
 AliFemtoCoulomb.cxx:245
 AliFemtoCoulomb.cxx:246
 AliFemtoCoulomb.cxx:247
 AliFemtoCoulomb.cxx:248
 AliFemtoCoulomb.cxx:249
 AliFemtoCoulomb.cxx:250
 AliFemtoCoulomb.cxx:251
 AliFemtoCoulomb.cxx:252
 AliFemtoCoulomb.cxx:253
 AliFemtoCoulomb.cxx:254
 AliFemtoCoulomb.cxx:255
 AliFemtoCoulomb.cxx:256
 AliFemtoCoulomb.cxx:257
 AliFemtoCoulomb.cxx:258
 AliFemtoCoulomb.cxx:259
 AliFemtoCoulomb.cxx:260
 AliFemtoCoulomb.cxx:261
 AliFemtoCoulomb.cxx:262
 AliFemtoCoulomb.cxx:263
 AliFemtoCoulomb.cxx:264
 AliFemtoCoulomb.cxx:265
 AliFemtoCoulomb.cxx:266
 AliFemtoCoulomb.cxx:267
 AliFemtoCoulomb.cxx:268
 AliFemtoCoulomb.cxx:269
 AliFemtoCoulomb.cxx:270
 AliFemtoCoulomb.cxx:271
 AliFemtoCoulomb.cxx:272
 AliFemtoCoulomb.cxx:273
 AliFemtoCoulomb.cxx:274
 AliFemtoCoulomb.cxx:275
 AliFemtoCoulomb.cxx:276
 AliFemtoCoulomb.cxx:277
 AliFemtoCoulomb.cxx:278
 AliFemtoCoulomb.cxx:279
 AliFemtoCoulomb.cxx:280
 AliFemtoCoulomb.cxx:281
 AliFemtoCoulomb.cxx:282
 AliFemtoCoulomb.cxx:283
 AliFemtoCoulomb.cxx:284
 AliFemtoCoulomb.cxx:285
 AliFemtoCoulomb.cxx:286
 AliFemtoCoulomb.cxx:287
 AliFemtoCoulomb.cxx:288
 AliFemtoCoulomb.cxx:289
 AliFemtoCoulomb.cxx:290
 AliFemtoCoulomb.cxx:291
 AliFemtoCoulomb.cxx:292
 AliFemtoCoulomb.cxx:293
 AliFemtoCoulomb.cxx:294
 AliFemtoCoulomb.cxx:295
 AliFemtoCoulomb.cxx:296
 AliFemtoCoulomb.cxx:297
 AliFemtoCoulomb.cxx:298
 AliFemtoCoulomb.cxx:299
 AliFemtoCoulomb.cxx:300
 AliFemtoCoulomb.cxx:301
 AliFemtoCoulomb.cxx:302
 AliFemtoCoulomb.cxx:303
 AliFemtoCoulomb.cxx:304
 AliFemtoCoulomb.cxx:305
 AliFemtoCoulomb.cxx:306
 AliFemtoCoulomb.cxx:307
 AliFemtoCoulomb.cxx:308
 AliFemtoCoulomb.cxx:309
 AliFemtoCoulomb.cxx:310
 AliFemtoCoulomb.cxx:311
 AliFemtoCoulomb.cxx:312
 AliFemtoCoulomb.cxx:313
 AliFemtoCoulomb.cxx:314
 AliFemtoCoulomb.cxx:315
 AliFemtoCoulomb.cxx:316
 AliFemtoCoulomb.cxx:317
 AliFemtoCoulomb.cxx:318
 AliFemtoCoulomb.cxx:319
 AliFemtoCoulomb.cxx:320
 AliFemtoCoulomb.cxx:321
 AliFemtoCoulomb.cxx:322
 AliFemtoCoulomb.cxx:323
 AliFemtoCoulomb.cxx:324
 AliFemtoCoulomb.cxx:325
 AliFemtoCoulomb.cxx:326
 AliFemtoCoulomb.cxx:327
 AliFemtoCoulomb.cxx:328
 AliFemtoCoulomb.cxx:329
 AliFemtoCoulomb.cxx:330
 AliFemtoCoulomb.cxx:331
 AliFemtoCoulomb.cxx:332
 AliFemtoCoulomb.cxx:333
 AliFemtoCoulomb.cxx:334
 AliFemtoCoulomb.cxx:335
 AliFemtoCoulomb.cxx:336
 AliFemtoCoulomb.cxx:337
 AliFemtoCoulomb.cxx:338
 AliFemtoCoulomb.cxx:339
 AliFemtoCoulomb.cxx:340
 AliFemtoCoulomb.cxx:341
 AliFemtoCoulomb.cxx:342
 AliFemtoCoulomb.cxx:343
 AliFemtoCoulomb.cxx:344
 AliFemtoCoulomb.cxx:345
 AliFemtoCoulomb.cxx:346
 AliFemtoCoulomb.cxx:347
 AliFemtoCoulomb.cxx:348
 AliFemtoCoulomb.cxx:349
 AliFemtoCoulomb.cxx:350
 AliFemtoCoulomb.cxx:351
 AliFemtoCoulomb.cxx:352
 AliFemtoCoulomb.cxx:353
 AliFemtoCoulomb.cxx:354
 AliFemtoCoulomb.cxx:355
 AliFemtoCoulomb.cxx:356
 AliFemtoCoulomb.cxx:357
 AliFemtoCoulomb.cxx:358
 AliFemtoCoulomb.cxx:359
 AliFemtoCoulomb.cxx:360
 AliFemtoCoulomb.cxx:361
 AliFemtoCoulomb.cxx:362
 AliFemtoCoulomb.cxx:363
 AliFemtoCoulomb.cxx:364
 AliFemtoCoulomb.cxx:365
 AliFemtoCoulomb.cxx:366
 AliFemtoCoulomb.cxx:367
 AliFemtoCoulomb.cxx:368
 AliFemtoCoulomb.cxx:369
 AliFemtoCoulomb.cxx:370
 AliFemtoCoulomb.cxx:371
 AliFemtoCoulomb.cxx:372
 AliFemtoCoulomb.cxx:373
 AliFemtoCoulomb.cxx:374
 AliFemtoCoulomb.cxx:375
 AliFemtoCoulomb.cxx:376
 AliFemtoCoulomb.cxx:377
 AliFemtoCoulomb.cxx:378
 AliFemtoCoulomb.cxx:379
 AliFemtoCoulomb.cxx:380
 AliFemtoCoulomb.cxx:381
 AliFemtoCoulomb.cxx:382
 AliFemtoCoulomb.cxx:383
 AliFemtoCoulomb.cxx:384
 AliFemtoCoulomb.cxx:385
 AliFemtoCoulomb.cxx:386
 AliFemtoCoulomb.cxx:387
 AliFemtoCoulomb.cxx:388
 AliFemtoCoulomb.cxx:389
 AliFemtoCoulomb.cxx:390
 AliFemtoCoulomb.cxx:391
 AliFemtoCoulomb.cxx:392
 AliFemtoCoulomb.cxx:393
 AliFemtoCoulomb.cxx:394
 AliFemtoCoulomb.cxx:395
 AliFemtoCoulomb.cxx:396
 AliFemtoCoulomb.cxx:397
 AliFemtoCoulomb.cxx:398
 AliFemtoCoulomb.cxx:399
 AliFemtoCoulomb.cxx:400
 AliFemtoCoulomb.cxx:401
 AliFemtoCoulomb.cxx:402
 AliFemtoCoulomb.cxx:403
 AliFemtoCoulomb.cxx:404
 AliFemtoCoulomb.cxx:405
 AliFemtoCoulomb.cxx:406
 AliFemtoCoulomb.cxx:407
 AliFemtoCoulomb.cxx:408
 AliFemtoCoulomb.cxx:409
 AliFemtoCoulomb.cxx:410
 AliFemtoCoulomb.cxx:411
 AliFemtoCoulomb.cxx:412
 AliFemtoCoulomb.cxx:413
 AliFemtoCoulomb.cxx:414
 AliFemtoCoulomb.cxx:415
 AliFemtoCoulomb.cxx:416
 AliFemtoCoulomb.cxx:417
 AliFemtoCoulomb.cxx:418
 AliFemtoCoulomb.cxx:419
 AliFemtoCoulomb.cxx:420
 AliFemtoCoulomb.cxx:421
 AliFemtoCoulomb.cxx:422
 AliFemtoCoulomb.cxx:423
 AliFemtoCoulomb.cxx:424
 AliFemtoCoulomb.cxx:425
 AliFemtoCoulomb.cxx:426
 AliFemtoCoulomb.cxx:427
 AliFemtoCoulomb.cxx:428
 AliFemtoCoulomb.cxx:429
 AliFemtoCoulomb.cxx:430
 AliFemtoCoulomb.cxx:431
 AliFemtoCoulomb.cxx:432
 AliFemtoCoulomb.cxx:433
 AliFemtoCoulomb.cxx:434
 AliFemtoCoulomb.cxx:435
 AliFemtoCoulomb.cxx:436
 AliFemtoCoulomb.cxx:437
 AliFemtoCoulomb.cxx:438
 AliFemtoCoulomb.cxx:439
 AliFemtoCoulomb.cxx:440
 AliFemtoCoulomb.cxx:441
 AliFemtoCoulomb.cxx:442
 AliFemtoCoulomb.cxx:443
 AliFemtoCoulomb.cxx:444
 AliFemtoCoulomb.cxx:445
 AliFemtoCoulomb.cxx:446
 AliFemtoCoulomb.cxx:447
 AliFemtoCoulomb.cxx:448
 AliFemtoCoulomb.cxx:449
 AliFemtoCoulomb.cxx:450
 AliFemtoCoulomb.cxx:451
 AliFemtoCoulomb.cxx:452
 AliFemtoCoulomb.cxx:453
 AliFemtoCoulomb.cxx:454
 AliFemtoCoulomb.cxx:455
 AliFemtoCoulomb.cxx:456
 AliFemtoCoulomb.cxx:457
 AliFemtoCoulomb.cxx:458
 AliFemtoCoulomb.cxx:459
 AliFemtoCoulomb.cxx:460
 AliFemtoCoulomb.cxx:461
 AliFemtoCoulomb.cxx:462
 AliFemtoCoulomb.cxx:463
 AliFemtoCoulomb.cxx:464
 AliFemtoCoulomb.cxx:465
 AliFemtoCoulomb.cxx:466
 AliFemtoCoulomb.cxx:467
 AliFemtoCoulomb.cxx:468
 AliFemtoCoulomb.cxx:469
 AliFemtoCoulomb.cxx:470
 AliFemtoCoulomb.cxx:471
 AliFemtoCoulomb.cxx:472