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

//
// AliGenGeVSim is a class implementing GeVSim event generator.
// 
// GeVSim is a simple Monte-Carlo event generator for testing detector and 
// algorythm performance especialy concerning flow and event-by-event studies
//
// In this event generator particles are generated from thermal distributions 
// without any dynamics and addicional constrains. Distribution parameters like
// multiplicity, particle type yields, inverse slope parameters, flow coeficients 
// and expansion velocities are expleicite defined by the user.
//
// GeVSim contains four thermal distributions the same as
// MevSim event generator developed for STAR experiment.
//
// In addition custom distributions can be used be the mean 
// either two dimensional formula (TF2), a two dimensional histogram or
// two one dimensional histograms.
//  
// Azimuthal distribution is deconvoluted from (Pt,Y) distribution
// and is described by two Fourier coefficients representing 
// Directed and Elliptic flow. 
// 
////////////////////////////////////////////////////////////////////////////////
//
// To apply flow to event ganerated by an arbitraly event generator
// refer to AliGenAfterBurnerFlow class.
//
////////////////////////////////////////////////////////////////////////////////
//
// For examples, parameters and testing macros refer to:
// http:/home.cern.ch/radomski
// 
// for more detailed description refer to ALICE NOTE
// "GeVSim Monte-Carlo Event Generator"
// S.Radosmki, P. Foka.
//  
// Author:
// Sylwester Radomski,
// GSI, March 2002
//  
// S.Radomski@gsi.de
//
////////////////////////////////////////////////////////////////////////////////
//
// Updated and revised: September 2002, S. Radomski, GSI
//
////////////////////////////////////////////////////////////////////////////////


#include <Riostream.h>
#include <TCanvas.h>
#include <TF1.h>
#include <TF2.h>
#include <TH1.h>
#include <TH2.h>
#include <TObjArray.h>
#include <TPDGCode.h>
#include <TParticle.h>
#include <TDatabasePDG.h>
#include <TROOT.h>


#include "AliGeVSimParticle.h"
#include "AliGenGeVSim.h"
#include "AliGenGeVSimEventHeader.h"
#include "AliGenerator.h"
#include "AliRun.h"


using std::cout;
using std::endl;
ClassImp(AliGenGeVSim)

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

AliGenGeVSim::AliGenGeVSim() : 
  AliGenerator(-1),
  fModel(0),
  fPsi(0),
  fIsMultTotal(kTRUE),
  fPtFormula(0),
  fYFormula(0),
  fPhiFormula(0),
  fCurrentForm(0),
  fPtYHist(0),
  fPartTypes(0) 
{
  //
  //  Default constructor
  // 

  for (Int_t i=0; i<4; i++)  
    fPtYFormula[i] = 0;
  for (Int_t i=0; i<2; i++)
    fHist[i] = 0;
}

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

AliGenGeVSim::AliGenGeVSim(Float_t psi, Bool_t isMultTotal) 
    : AliGenerator(-1),
      fModel(0),
      fPsi(psi),
      fIsMultTotal(isMultTotal),
      fPtFormula(0),
      fYFormula(0),
      fPhiFormula(0),
      fCurrentForm(0),
      fPtYHist(0),
      fPartTypes(0) 
 {
  //
  //  Standard Constructor.
  //  
  //  models - thermal model to be used:
  //          1    - deconvoluted pt and Y source
  //          2,3  - thermalized sphericaly symetric sources  
  //          4    - thermalized source with expansion
  //          5    - custom model defined in TF2 object named "gevsimPtY" 
  //          6    - custom model defined by two 1D histograms 
  //          7    - custom model defined by 2D histogram
  //
  //  psi   - reaction plane in degrees
  //  isMultTotal - multiplicity mode
  //                kTRUE - total multiplicity (default)
  //                kFALSE - dN/dY at midrapidity
  // 

  // checking consistancy
  
  if (psi < 0 || psi > 360 ) 
    Error ("AliGenGeVSim", "Reaction plane angle ( %13.3f )out of range [0..360]", psi);

  fPsi = psi * TMath::Pi() / 180. ;
  fIsMultTotal = isMultTotal;

  // Initialization 

  fPartTypes = new TObjArray();
  InitFormula();
}

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

AliGenGeVSim::~AliGenGeVSim() {
  //
  //  Default Destructor
  //  
  //  Removes TObjArray keeping list of registered particle types
  //

  if (fPartTypes != NULL) delete fPartTypes;
}


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

Bool_t AliGenGeVSim::CheckPtYPhi(Float_t pt, Float_t y, Float_t phi) const {
  //
  // private function used by Generate()
  //
  // Check bounds of Pt, Rapidity and Azimuthal Angle of a track
  // Used only when generating particles from a histogram
  //

  if ( TestBit(kPtRange) && ( pt < fPtMin || pt > fPtMax )) return kFALSE;
  if ( TestBit(kPhiRange) && ( phi < fPhiMin || phi > fPhiMax )) return kFALSE;
  if ( TestBit(kYRange) && ( y < fYMin || y > fYMax )) return kFALSE;

  return kTRUE;
}

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

Bool_t AliGenGeVSim::CheckAcceptance(Float_t p[3]) {
  //
  // private function used by Generate()
  //
  // Check bounds of a total momentum and theta of a track
  //

  if ( TestBit(kThetaRange) ) {
    
    Double_t theta = TMath::ATan2( TMath::Sqrt(p[0]*p[0]+p[1]*p[1]), p[2]);
    if ( theta < fThetaMin || theta > fThetaMax ) return kFALSE;
  }


  if ( TestBit(kMomentumRange) ) {
    
    Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
    if ( momentum < fPMin || momentum > fPMax) return kFALSE;
  }

  return kTRUE;
}

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

// Deconvoluted Pt Y formula

static Double_t aPtForm(Double_t * x, Double_t * par) {
  // ptForm: pt -> x[0] ,  mass -> [0] , temperature -> [1]
  // Description as string: " x * exp( -sqrt([0]*[0] + x*x) / [1] )"

    return x[0] * TMath::Exp( -sqrt(par[0]*par[0] + x[0]*x[0]) / par[1]);
  }

static  Double_t aYForm(Double_t * x, Double_t * par) {
  // y Form: y -> x[0] , sigmaY -> [0]
  // Description as string: " exp ( - x*x / (2 * [0]*[0] ) )"

    return TMath::Exp ( - x[0]*x[0] / (2 * par[0]*par[0] ) );
  }

// Models 1-3
// Description as strings:

//  const char *kFormE = " ( sqrt([0]*[0] + x*x) * cosh(y) ) ";
//  const char *kFormG = " ( 1 / sqrt( 1 - [2]*[2] ) ) ";
//  const char *kFormYp = "( [2]*sqrt(([0]*[0]+x*x)*cosh(y)*cosh(y)-[0]*[0])/([1]*sqrt(1-[2]*[2]))) ";

//  const char* kFormula[3] = {
//    " x * %s * exp( -%s / [1]) ", 
//    " (x * %s) / ( exp( %s / [1]) - 1 ) ",
//    " x*%s*exp(-%s*%s/[1])*((sinh(%s)/%s)+([1]/(%s*%s))*(sinh(%s)/%s-cosh(%s)))"
//  };
//   printf(kFormula[0], kFormE, kFormE);
//   printf(kFormula[1], kFormE, kFormE);
//   printf(kFormula[2], kFormE, kFormG, kFormE, kFormYp, kFormYp, kFormG, kFormE, kFormYp, kFormYp, kFormYp);


static Double_t aPtYFormula0(Double_t *x, Double_t * par) {
  // pt -> x , Y -> y
  // mass -> [0] , temperature -> [1] , expansion velocity -> [2]

  Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
  return x[0] * aFormE * TMath::Exp(-aFormE/par[1]);
}

static Double_t aPtYFormula1(Double_t *x, Double_t * par) {
  // pt -> x , Y -> y
  // mass -> [0] , temperature -> [1] , expansion velocity -> [2]

  Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
  return x[0] * aFormE / ( TMath::Exp( aFormE / par[1]) - 1 );
}

static Double_t aPtYFormula2(Double_t *x, Double_t * par) {
  // pt -> x , Y -> y
  // mass -> [0] , temperature -> [1] , expansion velocity -> [2]

  Double_t aFormE = TMath::Sqrt(par[0]*par[0] + x[0]*x[0]) * TMath::CosH(x[1]);
  Double_t aFormG = 1 / TMath::Sqrt((1.-par[2])*(1.+par[2]));
  Double_t aFormYp = par[2]*TMath::Sqrt( (par[0]*par[0] + x[0]*x[0]) 
					 * (TMath::CosH(x[1])-par[0])*(TMath::CosH(x[1])+par[0]))
    /( par[1]*TMath::Sqrt((1.-par[2])*(1.+par[2])));

  return x[0] * aFormE * TMath::Exp( - aFormG * aFormE / par[1])
    *( TMath::SinH(aFormYp)/aFormYp 
       + par[1]/(aFormG*aFormE) 
       * ( TMath::SinH(aFormYp)/aFormYp-TMath::CosH(aFormYp) ) );
}

// Phi Flow Formula

static Double_t aPhiForm(Double_t * x, Double_t * par) {
  // phi -> x
  // Psi -> [0] , Direct Flow -> [1] , Elliptical Flow -> [2]
  // Description as string: " 1 + 2*[1]*cos(x-[0]) + 2*[2]*cos(2*(x-[0])) "

  return 1 + 2*par[1]*TMath::Cos(x[0]-par[0]) 
    + 2*par[2]*TMath::Cos(2*(x[0]-par[0]));
}

void AliGenGeVSim::InitFormula() {
  //
  // private function
  //
  // Initalizes formulas used in GeVSim.

  // Deconvoluted Pt Y formula

  fPtFormula  = new TF1("gevsimPt", &aPtForm, 0, 3, 2);
  fYFormula   = new TF1("gevsimRapidity", &aYForm, -3, 3,1);

  fPtFormula->SetParNames("mass", "temperature");
  fPtFormula->SetParameters(1., 1.);
  
  fYFormula->SetParName(0, "sigmaY");
  fYFormula->SetParameter(0, 1.);

  // Grid for Pt and Y
  fPtFormula->SetNpx(100);
  fYFormula->SetNpx(100);
  

  // Models 1-3

  fPtYFormula[0] = new TF2("gevsimPtY_2", &aPtYFormula0, 0, 3, -2, 2, 2);

  fPtYFormula[1] = new TF2("gevsimPtY_3", &aPtYFormula1, 0, 3, -2, 2, 2);

  fPtYFormula[2] = new TF2("gevsimPtY_4", &aPtYFormula2, 0, 3, -2, 2, 3);

  fPtYFormula[3] = 0;


  // setting names & initialisation

  const char* kParamNames[3] = {"mass", "temperature", "expVel"};

  Int_t i, j;
  for (i=0; i<3; i++) {    

    fPtYFormula[i]->SetNpx(100);        // step 30 MeV  
    fPtYFormula[i]->SetNpy(100);        //

    for (j=0; j<3; j++) {

      if ( i != 2 && j == 2 ) continue; // ExpVel
      fPtYFormula[i]->SetParName(j, kParamNames[j]);
      fPtYFormula[i]->SetParameter(j, 0.5);
    }
  }
  
  // Phi Flow Formula

  fPhiFormula = new TF1("gevsimPhi", &aPhiForm, 0, 2*TMath::Pi(), 3);

  fPhiFormula->SetParNames("psi", "directed", "elliptic");
  fPhiFormula->SetParameters(0., 0., 0.);

  fPhiFormula->SetNpx(180);

}

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

void AliGenGeVSim::AddParticleType(AliGeVSimParticle *part) {
  //
  // Adds new type of particles.
  // 
  // Parameters are defeined in AliGeVSimParticle object
  // This method has to be called for every particle type
  //

  if (fPartTypes == NULL) 
     fPartTypes = new TObjArray();

  fPartTypes->Add(part);
}

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

void AliGenGeVSim::SetMultTotal(Bool_t isTotal) {
  //
  //
  //
  
  fIsMultTotal = isTotal;
}

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

Float_t AliGenGeVSim::FindScaler(Int_t paramId, Int_t pdg) {
  //
  // private function
  // Finds Scallar for a given parameter.
  // Function used in event-by-event mode.
  //
  // There are two types of scallars: deterministic and random
  // and they can work on either global or particle type level.
  // For every variable there are four scallars defined.  
  //  
  // Scallars are named as folowa
  // deterministic global level : "gevsimParam"        (eg. "gevsimTemp")
  // deterinistig type level    : "gevsimPdgParam"     (eg. "gevsim211Mult")
  // random global level        : "gevsimParamRndm"    (eg. "gevsimMultRndm")
  // random type level          : "gevsimPdgParamRndm" (eg. "gevsim-211V2Rndm");
  //
  // Pdg - code of a particle type in PDG standard (see: http://pdg.lbl.gov)
  // Param - parameter name. Allowed parameters:
  //
  // "Temp"   - inverse slope parameter
  // "SigmaY" - rapidity width - for model (1) only
  // "ExpVel" - expansion velocity - for model (4) only
  // "V1"     - directed flow
  // "V2"     - elliptic flow
  // "Mult"   - multiplicity
  //
  
  
  static const char* params[] = {"Temp", "SigmaY", "ExpVel", "V1", "V2", "Mult"};
  static const char* ending[] = {"", "Rndm"};

  static const char* patt1 = "gevsim%s%s";
  static const char* patt2 = "gevsim%d%s%s";

  char buffer[80];
  TF1 *form;
  
  Float_t scaler = 1.;

  // Scaler evoluation: i - global/local, j - determ/random

  Int_t i, j;

  for (i=0; i<2; i++) {
    for (j=0; j<2; j++) {
      
      form = 0;
      
      if (i == 0) snprintf(buffer, 80, patt1, params[paramId], ending[j]);      
      else snprintf(buffer, 80, patt2, pdg, params[paramId], ending[j]);
      
      form = (TF1 *)gROOT->GetFunction(buffer);

      if (form != 0) {
	if (j == 0) scaler *= form->Eval(gAlice->GetEvNumber()); 
	if (j == 1) {
	  form->SetParameter(0, gAlice->GetEvNumber());
	  scaler *= form->GetRandom();
	}
      }
    }
  }
  
  return scaler;
}

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

void AliGenGeVSim::DetermineReactionPlane() {
  //
  // private function used by Generate()
  //
  // Retermines Reaction Plane angle and set this value 
  // as a parameter [0] in fPhiFormula
  //
  // Note: if "gevsimPsiRndm" function is found it override both 
  //       "gevsimPhi" function and initial fPsi value
  //
  
  TF1 *form;
  
  form = 0;
  form = (TF1 *)gROOT->GetFunction("gevsimPsi");
  if (form) fPsi = form->Eval(gAlice->GetEvNumber()) * TMath::Pi() / 180;
  
  form = 0;
  form = (TF1 *)gROOT->GetFunction("gevsimPsiRndm");
  if (form) fPsi = form->GetRandom() * TMath::Pi() / 180;

  
  cout << "Psi = " << fPsi << "\t" << (Int_t)(fPsi*180./TMath::Pi()) << endl;
  
  fPhiFormula->SetParameter(0, fPsi);
}

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

Float_t AliGenGeVSim::GetdNdYToTotal() {
  //
  // Private, helper function used by Generate()
  //
  // Returns total multiplicity to dN/dY ration using current distribution.
  // The function have to be called after setting distribution and its 
  // parameters (like temperature).
  //

  Float_t integ, mag;
  const Double_t kMaxPt = 3.0, kMaxY = 2.; 

  if (fModel == 1) {
    
    integ = fYFormula->Integral(-kMaxY, kMaxY);
    mag = fYFormula->Eval(0);
    return integ/mag;
  }

  // 2D formula standard or custom

  if (fModel > 1 && fModel < 6) {
    
    integ =  ((TF2*)fCurrentForm)->Integral(0,kMaxPt, -kMaxY, kMaxY);
    mag = ((TF2*)fCurrentForm)->Integral(0, kMaxPt, -0.1, 0.1) / 0.2;
    return integ/mag;
  }

  // 2 1D histograms

  if (fModel == 6) {

    integ = fHist[1]->Integral(); 
    mag = fHist[0]->GetBinContent(fHist[0]->FindBin(0.));
    mag /= fHist[0]->GetBinWidth(fHist[0]->FindBin(0.));
    return integ/mag;
  }

  // 2D histogram
  
  if (fModel == 7) {

    // Not tested ...
    Int_t yBins = fPtYHist->GetNbinsY();
    Int_t ptBins = fPtYHist->GetNbinsX();

    integ = fPtYHist->Integral(0, ptBins, 0, yBins);
    mag = fPtYHist->Integral(0, ptBins, (yBins/2)-1, (yBins/2)+1 ) / 2;
    return integ/mag;
  }

  return 1;
}

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

void AliGenGeVSim::SetFormula(Int_t pdg) {
  //
  // Private function used by Generate() 
  //
  // Configure a formula for a given particle type and model Id (in fModel).
  // If custom formula or histogram was selected the function tries
  // to find it.
  //  
  // The function implements naming conventions for custom distributions names 
  // 

  char buff[40];
  const char* msg[4] = {
    "Custom Formula for Pt Y distribution not found [pdg = %d]",
    "Histogram for Pt distribution not found [pdg = %d]", 
    "Histogram for Y distribution not found [pdg = %d]",
    "HIstogram for Pt Y dostribution not found [pdg = %d]"
  };

  const char* pattern[8] = {
    "gevsimDistPtY", "gevsimDist%dPtY",
    "gevsimHistPt", "gevsimHist%dPt",
    "gevsimHistY", "gevsimHist%dY",
    "gevsimHistPtY", "gevsimHist%dPtY"
  };

  const char *where = "SetFormula";

  
  if (fModel < 1 || fModel > 7)
    Error("SetFormula", "Model Id (%d) out of range [1-7]", fModel);


  // standard models

  if (fModel == 1) fCurrentForm = fPtFormula;
  if (fModel > 1 && fModel < 5) fCurrentForm = fPtYFormula[fModel-2];


  // custom model defined by a formula 

  if (fModel == 5) {
    
    fCurrentForm = 0;
    fCurrentForm = (TF2*)gROOT->GetFunction(pattern[0]);
    
    if (!fCurrentForm) {

      snprintf(buff, 40, pattern[1], pdg);
      fCurrentForm = (TF2*)gROOT->GetFunction(buff);

      if (!fCurrentForm) Error(where, msg[0], pdg);
    }
  }

  // 2 1D histograms

  if (fModel == 6) {
    
    for (Int_t i=0; i<2; i++) {

      fHist[i] = 0;
      fHist[i] = (TH1D*)gROOT->FindObject(pattern[2+2*i]);
      
      if (!fHist[i]) {
	
	snprintf(buff, 40, pattern[3+2*i], pdg);
	fHist[i] = (TH1D*)gROOT->FindObject(buff);
	
	if (!fHist[i]) Error(where, msg[1+i], pdg);
      }
    }
  }
 
  // 2d histogram

  if (fModel == 7) {

    fPtYHist = 0;
    fPtYHist = (TH2D*)gROOT->FindObject(pattern[6]);
    
    if (!fPtYHist) {
      
      snprintf(buff, 40, pattern[7], pdg);
      fPtYHist = (TH2D*)gROOT->FindObject(buff);
    }

    if (!fPtYHist) Error(where, msg[3], pdg);
  }

}

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

void AliGenGeVSim:: AdjustFormula() {
  //
  // Private Function
  // Adjust fomula bounds according to acceptance cuts.
  //
  // Since GeVSim is producing "thermal" particles Pt
  // is cut at 3 GeV even when acceptance extends to grater momenta.
  //
  // WARNING !
  // If custom formula was provided function preserves
  // original cuts.
  //

  const Double_t kMaxPt = 3.0;
  const Double_t kMaxY = 2.0;
  Double_t minPt, maxPt, minY, maxY;

  
  if (fModel > 4) return;

  // max Pt 
  if (TestBit(kPtRange) && fPtMax < kMaxPt ) maxPt = fPtMax;
  else maxPt = kMaxPt;

  // min Pt
  if (TestBit(kPtRange)) minPt = fPtMin;
  else minPt = 0;

  if (TestBit(kPtRange) && fPtMin > kMaxPt ) 
    Warning("Acceptance", "Minimum Pt (%3.2f GeV) greater that 3.0 GeV ", fPtMin);

  // Max Pt < Max P
  if (TestBit(kMomentumRange) && fPtMax < maxPt) maxPt = fPtMax;
 
  // max and min rapidity
  if (TestBit(kYRange)) {
    minY = fYMin;
    maxY = fYMax;
  } else {
    minY = -kMaxY;
    maxY = kMaxY;
  }
  
  // adjust formula

  if (fModel == 1) {
    fPtFormula->SetRange(fPtMin, maxPt);
    fYFormula->SetRange(fYMin, fYMax);
  }
  
  if (fModel > 1)
    ((TF2*)fCurrentForm)->SetRange(minPt, minY, maxPt, maxY);

  // azimuthal cut

  if (TestBit(kPhiRange)) 
    fPhiFormula->SetRange(fPhiMin, fPhiMax);

}

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

void AliGenGeVSim::GetRandomPtY(Double_t &pt, Double_t &y) {
  //
  // Private function used by Generate()
  //
  // Returns random values of Pt and Y corresponding to selected
  // distribution.
  //
  
  if (fModel == 1) {
    pt = fPtFormula->GetRandom();
    y = fYFormula->GetRandom();
    return;
  }

  if (fModel > 1 && fModel < 6) {
    ((TF2*)fCurrentForm)->GetRandom2(pt, y);
    return;
  }
  
  if (fModel == 6) {
    pt = fHist[0]->GetRandom();
    y = fHist[1]->GetRandom();
  }
  
  if (fModel == 7) {
    fPtYHist->GetRandom2(pt, y);
    return;
  }
}

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

void AliGenGeVSim::Init() {
  //
  // Standard AliGenerator initializer.
  // does nothing
  //
}

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

void AliGenGeVSim::Generate() {
  //
  // Standard AliGenerator function
  // This function do actual job and puts particles on stack.
  //

  PDG_t pdg;                    // particle type
  Float_t mass;                 // particle mass
  Float_t orgin[3] = {0,0,0};   // particle orgin [cm]
  Float_t polar[3] = {0,0,0};   // polarisation
  Float_t time = 0;             // time of creation

  Float_t multiplicity = 0;           
  Bool_t isMultTotal = kTRUE;

  Float_t paramScaler;
  Float_t directedScaller = 1., ellipticScaller = 1.;

  TLorentzVector *v = new TLorentzVector(0,0,0,0);

  const Int_t kParent = -1;
  Int_t id;  

  // vertexing 
  VertexInternal();
  orgin[0] = fVertex[0];
  orgin[1] = fVertex[1];
  orgin[2] = fVertex[2];
  time = fTime;

  // Particle params database

  TDatabasePDG *db = TDatabasePDG::Instance(); 
  TParticlePDG *type;
  AliGeVSimParticle *partType;

  Int_t nType, nParticle, nParam;
  const Int_t kNParams = 6;

  // reaction plane determination and model
  DetermineReactionPlane();

  // loop over particle types

  for (nType = 0; nType < fPartTypes->GetEntries(); nType++) {

    partType = (AliGeVSimParticle *)fPartTypes->At(nType);

    pdg = (PDG_t)partType->GetPdgCode();
    type = db->GetParticle(pdg);
    mass = type->Mass();

    fModel = partType->GetModel();
    SetFormula(pdg);
    fCurrentForm->SetParameter("mass", mass);


    // Evaluation of parameters - loop over parameters

    for (nParam = 0; nParam < kNParams; nParam++) {
      
      paramScaler = FindScaler(nParam, pdg);

      if (nParam == 0)
	fCurrentForm->SetParameter("temperature", paramScaler * partType->GetTemperature());

      if (nParam == 1 && fModel == 1) 
	fYFormula->SetParameter("sigmaY", paramScaler * partType->GetSigmaY());
     
      if (nParam == 2 && fModel == 4) {

	Double_t totalExpVal =  paramScaler * partType->GetExpansionVelocity();

	if (totalExpVal == 0.0) totalExpVal = 0.0001;
	if (totalExpVal == 1.0) totalExpVal = 9.9999;

	fCurrentForm->SetParameter("expVel", totalExpVal);
      }

      // flow
     
      if (nParam == 3) directedScaller = paramScaler;
      if (nParam == 4) ellipticScaller = paramScaler;
      
      // multiplicity
      
      if (nParam == 5) {

	if (partType->IsMultForced()) isMultTotal = partType->IsMultTotal();
	else isMultTotal = fIsMultTotal;

	multiplicity = paramScaler * partType->GetMultiplicity();
	multiplicity *= (isMultTotal)? 1 : GetdNdYToTotal();
      } 
    }

    // Flow defined on the particle type level (not parameterised)
    if (partType->IsFlowSimple()) {
      fPhiFormula->SetParameter(1, partType->GetDirectedFlow(0,0) * directedScaller);
      fPhiFormula->SetParameter(2, partType->GetEllipticFlow(0,0) * ellipticScaller);
    }

    AdjustFormula();


    Info("Generate","PDG = %d \t Mult = %d", pdg, (Int_t)multiplicity);

    // loop over particles
    
    nParticle = 0;
    while (nParticle < multiplicity) {

      Double_t pt, y, phi;       // momentum in [pt,y,phi]
      Float_t p[3] = {0,0,0};    // particle momentum

      GetRandomPtY(pt, y);

      // phi distribution configuration when differential flow defined
      // to be optimised in future release 

      if (!partType->IsFlowSimple()) {
	fPhiFormula->SetParameter(1, partType->GetDirectedFlow(pt,y) * directedScaller);
	fPhiFormula->SetParameter(2, partType->GetEllipticFlow(pt,y) * ellipticScaller);
      }

      phi = fPhiFormula->GetRandom(); 

      if (!isMultTotal) nParticle++;
      if (fModel > 4 && !CheckPtYPhi(pt,y,phi) ) continue;
      
      // coordinate transformation
      v->SetPtEtaPhiM(pt, y, phi, mass);

      p[0] = v->Px();
      p[1] = v->Py();
      p[2] = v->Pz();

      // momentum range test
      if ( !CheckAcceptance(p) ) continue;

      // putting particle on the stack

      PushTrack(fTrackIt, kParent, pdg, p, orgin, polar, time, kPPrimary, id, fTrackIt);     
      if (isMultTotal) nParticle++;
    }
  }

  // prepare and store header

  AliGenGeVSimEventHeader *header = new AliGenGeVSimEventHeader("GeVSim header");
  TArrayF eventVertex(3,orgin);

  header->SetPrimaryVertex(eventVertex);
  header->SetInteractionTime(time);
  header->SetEventPlane(fPsi);
  header->SetEllipticFlow(fPhiFormula->GetParameter(2));

  gAlice->SetGenEventHeader(header);

  delete v;
}

//////////////////////////////////////////////////////////////////////////////////
 AliGenGeVSim.cxx:1
 AliGenGeVSim.cxx:2
 AliGenGeVSim.cxx:3
 AliGenGeVSim.cxx:4
 AliGenGeVSim.cxx:5
 AliGenGeVSim.cxx:6
 AliGenGeVSim.cxx:7
 AliGenGeVSim.cxx:8
 AliGenGeVSim.cxx:9
 AliGenGeVSim.cxx:10
 AliGenGeVSim.cxx:11
 AliGenGeVSim.cxx:12
 AliGenGeVSim.cxx:13
 AliGenGeVSim.cxx:14
 AliGenGeVSim.cxx:15
 AliGenGeVSim.cxx:16
 AliGenGeVSim.cxx:17
 AliGenGeVSim.cxx:18
 AliGenGeVSim.cxx:19
 AliGenGeVSim.cxx:20
 AliGenGeVSim.cxx:21
 AliGenGeVSim.cxx:22
 AliGenGeVSim.cxx:23
 AliGenGeVSim.cxx:24
 AliGenGeVSim.cxx:25
 AliGenGeVSim.cxx:26
 AliGenGeVSim.cxx:27
 AliGenGeVSim.cxx:28
 AliGenGeVSim.cxx:29
 AliGenGeVSim.cxx:30
 AliGenGeVSim.cxx:31
 AliGenGeVSim.cxx:32
 AliGenGeVSim.cxx:33
 AliGenGeVSim.cxx:34
 AliGenGeVSim.cxx:35
 AliGenGeVSim.cxx:36
 AliGenGeVSim.cxx:37
 AliGenGeVSim.cxx:38
 AliGenGeVSim.cxx:39
 AliGenGeVSim.cxx:40
 AliGenGeVSim.cxx:41
 AliGenGeVSim.cxx:42
 AliGenGeVSim.cxx:43
 AliGenGeVSim.cxx:44
 AliGenGeVSim.cxx:45
 AliGenGeVSim.cxx:46
 AliGenGeVSim.cxx:47
 AliGenGeVSim.cxx:48
 AliGenGeVSim.cxx:49
 AliGenGeVSim.cxx:50
 AliGenGeVSim.cxx:51
 AliGenGeVSim.cxx:52
 AliGenGeVSim.cxx:53
 AliGenGeVSim.cxx:54
 AliGenGeVSim.cxx:55
 AliGenGeVSim.cxx:56
 AliGenGeVSim.cxx:57
 AliGenGeVSim.cxx:58
 AliGenGeVSim.cxx:59
 AliGenGeVSim.cxx:60
 AliGenGeVSim.cxx:61
 AliGenGeVSim.cxx:62
 AliGenGeVSim.cxx:63
 AliGenGeVSim.cxx:64
 AliGenGeVSim.cxx:65
 AliGenGeVSim.cxx:66
 AliGenGeVSim.cxx:67
 AliGenGeVSim.cxx:68
 AliGenGeVSim.cxx:69
 AliGenGeVSim.cxx:70
 AliGenGeVSim.cxx:71
 AliGenGeVSim.cxx:72
 AliGenGeVSim.cxx:73
 AliGenGeVSim.cxx:74
 AliGenGeVSim.cxx:75
 AliGenGeVSim.cxx:76
 AliGenGeVSim.cxx:77
 AliGenGeVSim.cxx:78
 AliGenGeVSim.cxx:79
 AliGenGeVSim.cxx:80
 AliGenGeVSim.cxx:81
 AliGenGeVSim.cxx:82
 AliGenGeVSim.cxx:83
 AliGenGeVSim.cxx:84
 AliGenGeVSim.cxx:85
 AliGenGeVSim.cxx:86
 AliGenGeVSim.cxx:87
 AliGenGeVSim.cxx:88
 AliGenGeVSim.cxx:89
 AliGenGeVSim.cxx:90
 AliGenGeVSim.cxx:91
 AliGenGeVSim.cxx:92
 AliGenGeVSim.cxx:93
 AliGenGeVSim.cxx:94
 AliGenGeVSim.cxx:95
 AliGenGeVSim.cxx:96
 AliGenGeVSim.cxx:97
 AliGenGeVSim.cxx:98
 AliGenGeVSim.cxx:99
 AliGenGeVSim.cxx:100
 AliGenGeVSim.cxx:101
 AliGenGeVSim.cxx:102
 AliGenGeVSim.cxx:103
 AliGenGeVSim.cxx:104
 AliGenGeVSim.cxx:105
 AliGenGeVSim.cxx:106
 AliGenGeVSim.cxx:107
 AliGenGeVSim.cxx:108
 AliGenGeVSim.cxx:109
 AliGenGeVSim.cxx:110
 AliGenGeVSim.cxx:111
 AliGenGeVSim.cxx:112
 AliGenGeVSim.cxx:113
 AliGenGeVSim.cxx:114
 AliGenGeVSim.cxx:115
 AliGenGeVSim.cxx:116
 AliGenGeVSim.cxx:117
 AliGenGeVSim.cxx:118
 AliGenGeVSim.cxx:119
 AliGenGeVSim.cxx:120
 AliGenGeVSim.cxx:121
 AliGenGeVSim.cxx:122
 AliGenGeVSim.cxx:123
 AliGenGeVSim.cxx:124
 AliGenGeVSim.cxx:125
 AliGenGeVSim.cxx:126
 AliGenGeVSim.cxx:127
 AliGenGeVSim.cxx:128
 AliGenGeVSim.cxx:129
 AliGenGeVSim.cxx:130
 AliGenGeVSim.cxx:131
 AliGenGeVSim.cxx:132
 AliGenGeVSim.cxx:133
 AliGenGeVSim.cxx:134
 AliGenGeVSim.cxx:135
 AliGenGeVSim.cxx:136
 AliGenGeVSim.cxx:137
 AliGenGeVSim.cxx:138
 AliGenGeVSim.cxx:139
 AliGenGeVSim.cxx:140
 AliGenGeVSim.cxx:141
 AliGenGeVSim.cxx:142
 AliGenGeVSim.cxx:143
 AliGenGeVSim.cxx:144
 AliGenGeVSim.cxx:145
 AliGenGeVSim.cxx:146
 AliGenGeVSim.cxx:147
 AliGenGeVSim.cxx:148
 AliGenGeVSim.cxx:149
 AliGenGeVSim.cxx:150
 AliGenGeVSim.cxx:151
 AliGenGeVSim.cxx:152
 AliGenGeVSim.cxx:153
 AliGenGeVSim.cxx:154
 AliGenGeVSim.cxx:155
 AliGenGeVSim.cxx:156
 AliGenGeVSim.cxx:157
 AliGenGeVSim.cxx:158
 AliGenGeVSim.cxx:159
 AliGenGeVSim.cxx:160
 AliGenGeVSim.cxx:161
 AliGenGeVSim.cxx:162
 AliGenGeVSim.cxx:163
 AliGenGeVSim.cxx:164
 AliGenGeVSim.cxx:165
 AliGenGeVSim.cxx:166
 AliGenGeVSim.cxx:167
 AliGenGeVSim.cxx:168
 AliGenGeVSim.cxx:169
 AliGenGeVSim.cxx:170
 AliGenGeVSim.cxx:171
 AliGenGeVSim.cxx:172
 AliGenGeVSim.cxx:173
 AliGenGeVSim.cxx:174
 AliGenGeVSim.cxx:175
 AliGenGeVSim.cxx:176
 AliGenGeVSim.cxx:177
 AliGenGeVSim.cxx:178
 AliGenGeVSim.cxx:179
 AliGenGeVSim.cxx:180
 AliGenGeVSim.cxx:181
 AliGenGeVSim.cxx:182
 AliGenGeVSim.cxx:183
 AliGenGeVSim.cxx:184
 AliGenGeVSim.cxx:185
 AliGenGeVSim.cxx:186
 AliGenGeVSim.cxx:187
 AliGenGeVSim.cxx:188
 AliGenGeVSim.cxx:189
 AliGenGeVSim.cxx:190
 AliGenGeVSim.cxx:191
 AliGenGeVSim.cxx:192
 AliGenGeVSim.cxx:193
 AliGenGeVSim.cxx:194
 AliGenGeVSim.cxx:195
 AliGenGeVSim.cxx:196
 AliGenGeVSim.cxx:197
 AliGenGeVSim.cxx:198
 AliGenGeVSim.cxx:199
 AliGenGeVSim.cxx:200
 AliGenGeVSim.cxx:201
 AliGenGeVSim.cxx:202
 AliGenGeVSim.cxx:203
 AliGenGeVSim.cxx:204
 AliGenGeVSim.cxx:205
 AliGenGeVSim.cxx:206
 AliGenGeVSim.cxx:207
 AliGenGeVSim.cxx:208
 AliGenGeVSim.cxx:209
 AliGenGeVSim.cxx:210
 AliGenGeVSim.cxx:211
 AliGenGeVSim.cxx:212
 AliGenGeVSim.cxx:213
 AliGenGeVSim.cxx:214
 AliGenGeVSim.cxx:215
 AliGenGeVSim.cxx:216
 AliGenGeVSim.cxx:217
 AliGenGeVSim.cxx:218
 AliGenGeVSim.cxx:219
 AliGenGeVSim.cxx:220
 AliGenGeVSim.cxx:221
 AliGenGeVSim.cxx:222
 AliGenGeVSim.cxx:223
 AliGenGeVSim.cxx:224
 AliGenGeVSim.cxx:225
 AliGenGeVSim.cxx:226
 AliGenGeVSim.cxx:227
 AliGenGeVSim.cxx:228
 AliGenGeVSim.cxx:229
 AliGenGeVSim.cxx:230
 AliGenGeVSim.cxx:231
 AliGenGeVSim.cxx:232
 AliGenGeVSim.cxx:233
 AliGenGeVSim.cxx:234
 AliGenGeVSim.cxx:235
 AliGenGeVSim.cxx:236
 AliGenGeVSim.cxx:237
 AliGenGeVSim.cxx:238
 AliGenGeVSim.cxx:239
 AliGenGeVSim.cxx:240
 AliGenGeVSim.cxx:241
 AliGenGeVSim.cxx:242
 AliGenGeVSim.cxx:243
 AliGenGeVSim.cxx:244
 AliGenGeVSim.cxx:245
 AliGenGeVSim.cxx:246
 AliGenGeVSim.cxx:247
 AliGenGeVSim.cxx:248
 AliGenGeVSim.cxx:249
 AliGenGeVSim.cxx:250
 AliGenGeVSim.cxx:251
 AliGenGeVSim.cxx:252
 AliGenGeVSim.cxx:253
 AliGenGeVSim.cxx:254
 AliGenGeVSim.cxx:255
 AliGenGeVSim.cxx:256
 AliGenGeVSim.cxx:257
 AliGenGeVSim.cxx:258
 AliGenGeVSim.cxx:259
 AliGenGeVSim.cxx:260
 AliGenGeVSim.cxx:261
 AliGenGeVSim.cxx:262
 AliGenGeVSim.cxx:263
 AliGenGeVSim.cxx:264
 AliGenGeVSim.cxx:265
 AliGenGeVSim.cxx:266
 AliGenGeVSim.cxx:267
 AliGenGeVSim.cxx:268
 AliGenGeVSim.cxx:269
 AliGenGeVSim.cxx:270
 AliGenGeVSim.cxx:271
 AliGenGeVSim.cxx:272
 AliGenGeVSim.cxx:273
 AliGenGeVSim.cxx:274
 AliGenGeVSim.cxx:275
 AliGenGeVSim.cxx:276
 AliGenGeVSim.cxx:277
 AliGenGeVSim.cxx:278
 AliGenGeVSim.cxx:279
 AliGenGeVSim.cxx:280
 AliGenGeVSim.cxx:281
 AliGenGeVSim.cxx:282
 AliGenGeVSim.cxx:283
 AliGenGeVSim.cxx:284
 AliGenGeVSim.cxx:285
 AliGenGeVSim.cxx:286
 AliGenGeVSim.cxx:287
 AliGenGeVSim.cxx:288
 AliGenGeVSim.cxx:289
 AliGenGeVSim.cxx:290
 AliGenGeVSim.cxx:291
 AliGenGeVSim.cxx:292
 AliGenGeVSim.cxx:293
 AliGenGeVSim.cxx:294
 AliGenGeVSim.cxx:295
 AliGenGeVSim.cxx:296
 AliGenGeVSim.cxx:297
 AliGenGeVSim.cxx:298
 AliGenGeVSim.cxx:299
 AliGenGeVSim.cxx:300
 AliGenGeVSim.cxx:301
 AliGenGeVSim.cxx:302
 AliGenGeVSim.cxx:303
 AliGenGeVSim.cxx:304
 AliGenGeVSim.cxx:305
 AliGenGeVSim.cxx:306
 AliGenGeVSim.cxx:307
 AliGenGeVSim.cxx:308
 AliGenGeVSim.cxx:309
 AliGenGeVSim.cxx:310
 AliGenGeVSim.cxx:311
 AliGenGeVSim.cxx:312
 AliGenGeVSim.cxx:313
 AliGenGeVSim.cxx:314
 AliGenGeVSim.cxx:315
 AliGenGeVSim.cxx:316
 AliGenGeVSim.cxx:317
 AliGenGeVSim.cxx:318
 AliGenGeVSim.cxx:319
 AliGenGeVSim.cxx:320
 AliGenGeVSim.cxx:321
 AliGenGeVSim.cxx:322
 AliGenGeVSim.cxx:323
 AliGenGeVSim.cxx:324
 AliGenGeVSim.cxx:325
 AliGenGeVSim.cxx:326
 AliGenGeVSim.cxx:327
 AliGenGeVSim.cxx:328
 AliGenGeVSim.cxx:329
 AliGenGeVSim.cxx:330
 AliGenGeVSim.cxx:331
 AliGenGeVSim.cxx:332
 AliGenGeVSim.cxx:333
 AliGenGeVSim.cxx:334
 AliGenGeVSim.cxx:335
 AliGenGeVSim.cxx:336
 AliGenGeVSim.cxx:337
 AliGenGeVSim.cxx:338
 AliGenGeVSim.cxx:339
 AliGenGeVSim.cxx:340
 AliGenGeVSim.cxx:341
 AliGenGeVSim.cxx:342
 AliGenGeVSim.cxx:343
 AliGenGeVSim.cxx:344
 AliGenGeVSim.cxx:345
 AliGenGeVSim.cxx:346
 AliGenGeVSim.cxx:347
 AliGenGeVSim.cxx:348
 AliGenGeVSim.cxx:349
 AliGenGeVSim.cxx:350
 AliGenGeVSim.cxx:351
 AliGenGeVSim.cxx:352
 AliGenGeVSim.cxx:353
 AliGenGeVSim.cxx:354
 AliGenGeVSim.cxx:355
 AliGenGeVSim.cxx:356
 AliGenGeVSim.cxx:357
 AliGenGeVSim.cxx:358
 AliGenGeVSim.cxx:359
 AliGenGeVSim.cxx:360
 AliGenGeVSim.cxx:361
 AliGenGeVSim.cxx:362
 AliGenGeVSim.cxx:363
 AliGenGeVSim.cxx:364
 AliGenGeVSim.cxx:365
 AliGenGeVSim.cxx:366
 AliGenGeVSim.cxx:367
 AliGenGeVSim.cxx:368
 AliGenGeVSim.cxx:369
 AliGenGeVSim.cxx:370
 AliGenGeVSim.cxx:371
 AliGenGeVSim.cxx:372
 AliGenGeVSim.cxx:373
 AliGenGeVSim.cxx:374
 AliGenGeVSim.cxx:375
 AliGenGeVSim.cxx:376
 AliGenGeVSim.cxx:377
 AliGenGeVSim.cxx:378
 AliGenGeVSim.cxx:379
 AliGenGeVSim.cxx:380
 AliGenGeVSim.cxx:381
 AliGenGeVSim.cxx:382
 AliGenGeVSim.cxx:383
 AliGenGeVSim.cxx:384
 AliGenGeVSim.cxx:385
 AliGenGeVSim.cxx:386
 AliGenGeVSim.cxx:387
 AliGenGeVSim.cxx:388
 AliGenGeVSim.cxx:389
 AliGenGeVSim.cxx:390
 AliGenGeVSim.cxx:391
 AliGenGeVSim.cxx:392
 AliGenGeVSim.cxx:393
 AliGenGeVSim.cxx:394
 AliGenGeVSim.cxx:395
 AliGenGeVSim.cxx:396
 AliGenGeVSim.cxx:397
 AliGenGeVSim.cxx:398
 AliGenGeVSim.cxx:399
 AliGenGeVSim.cxx:400
 AliGenGeVSim.cxx:401
 AliGenGeVSim.cxx:402
 AliGenGeVSim.cxx:403
 AliGenGeVSim.cxx:404
 AliGenGeVSim.cxx:405
 AliGenGeVSim.cxx:406
 AliGenGeVSim.cxx:407
 AliGenGeVSim.cxx:408
 AliGenGeVSim.cxx:409
 AliGenGeVSim.cxx:410
 AliGenGeVSim.cxx:411
 AliGenGeVSim.cxx:412
 AliGenGeVSim.cxx:413
 AliGenGeVSim.cxx:414
 AliGenGeVSim.cxx:415
 AliGenGeVSim.cxx:416
 AliGenGeVSim.cxx:417
 AliGenGeVSim.cxx:418
 AliGenGeVSim.cxx:419
 AliGenGeVSim.cxx:420
 AliGenGeVSim.cxx:421
 AliGenGeVSim.cxx:422
 AliGenGeVSim.cxx:423
 AliGenGeVSim.cxx:424
 AliGenGeVSim.cxx:425
 AliGenGeVSim.cxx:426
 AliGenGeVSim.cxx:427
 AliGenGeVSim.cxx:428
 AliGenGeVSim.cxx:429
 AliGenGeVSim.cxx:430
 AliGenGeVSim.cxx:431
 AliGenGeVSim.cxx:432
 AliGenGeVSim.cxx:433
 AliGenGeVSim.cxx:434
 AliGenGeVSim.cxx:435
 AliGenGeVSim.cxx:436
 AliGenGeVSim.cxx:437
 AliGenGeVSim.cxx:438
 AliGenGeVSim.cxx:439
 AliGenGeVSim.cxx:440
 AliGenGeVSim.cxx:441
 AliGenGeVSim.cxx:442
 AliGenGeVSim.cxx:443
 AliGenGeVSim.cxx:444
 AliGenGeVSim.cxx:445
 AliGenGeVSim.cxx:446
 AliGenGeVSim.cxx:447
 AliGenGeVSim.cxx:448
 AliGenGeVSim.cxx:449
 AliGenGeVSim.cxx:450
 AliGenGeVSim.cxx:451
 AliGenGeVSim.cxx:452
 AliGenGeVSim.cxx:453
 AliGenGeVSim.cxx:454
 AliGenGeVSim.cxx:455
 AliGenGeVSim.cxx:456
 AliGenGeVSim.cxx:457
 AliGenGeVSim.cxx:458
 AliGenGeVSim.cxx:459
 AliGenGeVSim.cxx:460
 AliGenGeVSim.cxx:461
 AliGenGeVSim.cxx:462
 AliGenGeVSim.cxx:463
 AliGenGeVSim.cxx:464
 AliGenGeVSim.cxx:465
 AliGenGeVSim.cxx:466
 AliGenGeVSim.cxx:467
 AliGenGeVSim.cxx:468
 AliGenGeVSim.cxx:469
 AliGenGeVSim.cxx:470
 AliGenGeVSim.cxx:471
 AliGenGeVSim.cxx:472
 AliGenGeVSim.cxx:473
 AliGenGeVSim.cxx:474
 AliGenGeVSim.cxx:475
 AliGenGeVSim.cxx:476
 AliGenGeVSim.cxx:477
 AliGenGeVSim.cxx:478
 AliGenGeVSim.cxx:479
 AliGenGeVSim.cxx:480
 AliGenGeVSim.cxx:481
 AliGenGeVSim.cxx:482
 AliGenGeVSim.cxx:483
 AliGenGeVSim.cxx:484
 AliGenGeVSim.cxx:485
 AliGenGeVSim.cxx:486
 AliGenGeVSim.cxx:487
 AliGenGeVSim.cxx:488
 AliGenGeVSim.cxx:489
 AliGenGeVSim.cxx:490
 AliGenGeVSim.cxx:491
 AliGenGeVSim.cxx:492
 AliGenGeVSim.cxx:493
 AliGenGeVSim.cxx:494
 AliGenGeVSim.cxx:495
 AliGenGeVSim.cxx:496
 AliGenGeVSim.cxx:497
 AliGenGeVSim.cxx:498
 AliGenGeVSim.cxx:499
 AliGenGeVSim.cxx:500
 AliGenGeVSim.cxx:501
 AliGenGeVSim.cxx:502
 AliGenGeVSim.cxx:503
 AliGenGeVSim.cxx:504
 AliGenGeVSim.cxx:505
 AliGenGeVSim.cxx:506
 AliGenGeVSim.cxx:507
 AliGenGeVSim.cxx:508
 AliGenGeVSim.cxx:509
 AliGenGeVSim.cxx:510
 AliGenGeVSim.cxx:511
 AliGenGeVSim.cxx:512
 AliGenGeVSim.cxx:513
 AliGenGeVSim.cxx:514
 AliGenGeVSim.cxx:515
 AliGenGeVSim.cxx:516
 AliGenGeVSim.cxx:517
 AliGenGeVSim.cxx:518
 AliGenGeVSim.cxx:519
 AliGenGeVSim.cxx:520
 AliGenGeVSim.cxx:521
 AliGenGeVSim.cxx:522
 AliGenGeVSim.cxx:523
 AliGenGeVSim.cxx:524
 AliGenGeVSim.cxx:525
 AliGenGeVSim.cxx:526
 AliGenGeVSim.cxx:527
 AliGenGeVSim.cxx:528
 AliGenGeVSim.cxx:529
 AliGenGeVSim.cxx:530
 AliGenGeVSim.cxx:531
 AliGenGeVSim.cxx:532
 AliGenGeVSim.cxx:533
 AliGenGeVSim.cxx:534
 AliGenGeVSim.cxx:535
 AliGenGeVSim.cxx:536
 AliGenGeVSim.cxx:537
 AliGenGeVSim.cxx:538
 AliGenGeVSim.cxx:539
 AliGenGeVSim.cxx:540
 AliGenGeVSim.cxx:541
 AliGenGeVSim.cxx:542
 AliGenGeVSim.cxx:543
 AliGenGeVSim.cxx:544
 AliGenGeVSim.cxx:545
 AliGenGeVSim.cxx:546
 AliGenGeVSim.cxx:547
 AliGenGeVSim.cxx:548
 AliGenGeVSim.cxx:549
 AliGenGeVSim.cxx:550
 AliGenGeVSim.cxx:551
 AliGenGeVSim.cxx:552
 AliGenGeVSim.cxx:553
 AliGenGeVSim.cxx:554
 AliGenGeVSim.cxx:555
 AliGenGeVSim.cxx:556
 AliGenGeVSim.cxx:557
 AliGenGeVSim.cxx:558
 AliGenGeVSim.cxx:559
 AliGenGeVSim.cxx:560
 AliGenGeVSim.cxx:561
 AliGenGeVSim.cxx:562
 AliGenGeVSim.cxx:563
 AliGenGeVSim.cxx:564
 AliGenGeVSim.cxx:565
 AliGenGeVSim.cxx:566
 AliGenGeVSim.cxx:567
 AliGenGeVSim.cxx:568
 AliGenGeVSim.cxx:569
 AliGenGeVSim.cxx:570
 AliGenGeVSim.cxx:571
 AliGenGeVSim.cxx:572
 AliGenGeVSim.cxx:573
 AliGenGeVSim.cxx:574
 AliGenGeVSim.cxx:575
 AliGenGeVSim.cxx:576
 AliGenGeVSim.cxx:577
 AliGenGeVSim.cxx:578
 AliGenGeVSim.cxx:579
 AliGenGeVSim.cxx:580
 AliGenGeVSim.cxx:581
 AliGenGeVSim.cxx:582
 AliGenGeVSim.cxx:583
 AliGenGeVSim.cxx:584
 AliGenGeVSim.cxx:585
 AliGenGeVSim.cxx:586
 AliGenGeVSim.cxx:587
 AliGenGeVSim.cxx:588
 AliGenGeVSim.cxx:589
 AliGenGeVSim.cxx:590
 AliGenGeVSim.cxx:591
 AliGenGeVSim.cxx:592
 AliGenGeVSim.cxx:593
 AliGenGeVSim.cxx:594
 AliGenGeVSim.cxx:595
 AliGenGeVSim.cxx:596
 AliGenGeVSim.cxx:597
 AliGenGeVSim.cxx:598
 AliGenGeVSim.cxx:599
 AliGenGeVSim.cxx:600
 AliGenGeVSim.cxx:601
 AliGenGeVSim.cxx:602
 AliGenGeVSim.cxx:603
 AliGenGeVSim.cxx:604
 AliGenGeVSim.cxx:605
 AliGenGeVSim.cxx:606
 AliGenGeVSim.cxx:607
 AliGenGeVSim.cxx:608
 AliGenGeVSim.cxx:609
 AliGenGeVSim.cxx:610
 AliGenGeVSim.cxx:611
 AliGenGeVSim.cxx:612
 AliGenGeVSim.cxx:613
 AliGenGeVSim.cxx:614
 AliGenGeVSim.cxx:615
 AliGenGeVSim.cxx:616
 AliGenGeVSim.cxx:617
 AliGenGeVSim.cxx:618
 AliGenGeVSim.cxx:619
 AliGenGeVSim.cxx:620
 AliGenGeVSim.cxx:621
 AliGenGeVSim.cxx:622
 AliGenGeVSim.cxx:623
 AliGenGeVSim.cxx:624
 AliGenGeVSim.cxx:625
 AliGenGeVSim.cxx:626
 AliGenGeVSim.cxx:627
 AliGenGeVSim.cxx:628
 AliGenGeVSim.cxx:629
 AliGenGeVSim.cxx:630
 AliGenGeVSim.cxx:631
 AliGenGeVSim.cxx:632
 AliGenGeVSim.cxx:633
 AliGenGeVSim.cxx:634
 AliGenGeVSim.cxx:635
 AliGenGeVSim.cxx:636
 AliGenGeVSim.cxx:637
 AliGenGeVSim.cxx:638
 AliGenGeVSim.cxx:639
 AliGenGeVSim.cxx:640
 AliGenGeVSim.cxx:641
 AliGenGeVSim.cxx:642
 AliGenGeVSim.cxx:643
 AliGenGeVSim.cxx:644
 AliGenGeVSim.cxx:645
 AliGenGeVSim.cxx:646
 AliGenGeVSim.cxx:647
 AliGenGeVSim.cxx:648
 AliGenGeVSim.cxx:649
 AliGenGeVSim.cxx:650
 AliGenGeVSim.cxx:651
 AliGenGeVSim.cxx:652
 AliGenGeVSim.cxx:653
 AliGenGeVSim.cxx:654
 AliGenGeVSim.cxx:655
 AliGenGeVSim.cxx:656
 AliGenGeVSim.cxx:657
 AliGenGeVSim.cxx:658
 AliGenGeVSim.cxx:659
 AliGenGeVSim.cxx:660
 AliGenGeVSim.cxx:661
 AliGenGeVSim.cxx:662
 AliGenGeVSim.cxx:663
 AliGenGeVSim.cxx:664
 AliGenGeVSim.cxx:665
 AliGenGeVSim.cxx:666
 AliGenGeVSim.cxx:667
 AliGenGeVSim.cxx:668
 AliGenGeVSim.cxx:669
 AliGenGeVSim.cxx:670
 AliGenGeVSim.cxx:671
 AliGenGeVSim.cxx:672
 AliGenGeVSim.cxx:673
 AliGenGeVSim.cxx:674
 AliGenGeVSim.cxx:675
 AliGenGeVSim.cxx:676
 AliGenGeVSim.cxx:677
 AliGenGeVSim.cxx:678
 AliGenGeVSim.cxx:679
 AliGenGeVSim.cxx:680
 AliGenGeVSim.cxx:681
 AliGenGeVSim.cxx:682
 AliGenGeVSim.cxx:683
 AliGenGeVSim.cxx:684
 AliGenGeVSim.cxx:685
 AliGenGeVSim.cxx:686
 AliGenGeVSim.cxx:687
 AliGenGeVSim.cxx:688
 AliGenGeVSim.cxx:689
 AliGenGeVSim.cxx:690
 AliGenGeVSim.cxx:691
 AliGenGeVSim.cxx:692
 AliGenGeVSim.cxx:693
 AliGenGeVSim.cxx:694
 AliGenGeVSim.cxx:695
 AliGenGeVSim.cxx:696
 AliGenGeVSim.cxx:697
 AliGenGeVSim.cxx:698
 AliGenGeVSim.cxx:699
 AliGenGeVSim.cxx:700
 AliGenGeVSim.cxx:701
 AliGenGeVSim.cxx:702
 AliGenGeVSim.cxx:703
 AliGenGeVSim.cxx:704
 AliGenGeVSim.cxx:705
 AliGenGeVSim.cxx:706
 AliGenGeVSim.cxx:707
 AliGenGeVSim.cxx:708
 AliGenGeVSim.cxx:709
 AliGenGeVSim.cxx:710
 AliGenGeVSim.cxx:711
 AliGenGeVSim.cxx:712
 AliGenGeVSim.cxx:713
 AliGenGeVSim.cxx:714
 AliGenGeVSim.cxx:715
 AliGenGeVSim.cxx:716
 AliGenGeVSim.cxx:717
 AliGenGeVSim.cxx:718
 AliGenGeVSim.cxx:719
 AliGenGeVSim.cxx:720
 AliGenGeVSim.cxx:721
 AliGenGeVSim.cxx:722
 AliGenGeVSim.cxx:723
 AliGenGeVSim.cxx:724
 AliGenGeVSim.cxx:725
 AliGenGeVSim.cxx:726
 AliGenGeVSim.cxx:727
 AliGenGeVSim.cxx:728
 AliGenGeVSim.cxx:729
 AliGenGeVSim.cxx:730
 AliGenGeVSim.cxx:731
 AliGenGeVSim.cxx:732
 AliGenGeVSim.cxx:733
 AliGenGeVSim.cxx:734
 AliGenGeVSim.cxx:735
 AliGenGeVSim.cxx:736
 AliGenGeVSim.cxx:737
 AliGenGeVSim.cxx:738
 AliGenGeVSim.cxx:739
 AliGenGeVSim.cxx:740
 AliGenGeVSim.cxx:741
 AliGenGeVSim.cxx:742
 AliGenGeVSim.cxx:743
 AliGenGeVSim.cxx:744
 AliGenGeVSim.cxx:745
 AliGenGeVSim.cxx:746
 AliGenGeVSim.cxx:747
 AliGenGeVSim.cxx:748
 AliGenGeVSim.cxx:749
 AliGenGeVSim.cxx:750
 AliGenGeVSim.cxx:751
 AliGenGeVSim.cxx:752
 AliGenGeVSim.cxx:753
 AliGenGeVSim.cxx:754
 AliGenGeVSim.cxx:755
 AliGenGeVSim.cxx:756
 AliGenGeVSim.cxx:757
 AliGenGeVSim.cxx:758
 AliGenGeVSim.cxx:759
 AliGenGeVSim.cxx:760
 AliGenGeVSim.cxx:761
 AliGenGeVSim.cxx:762
 AliGenGeVSim.cxx:763
 AliGenGeVSim.cxx:764
 AliGenGeVSim.cxx:765
 AliGenGeVSim.cxx:766
 AliGenGeVSim.cxx:767
 AliGenGeVSim.cxx:768
 AliGenGeVSim.cxx:769
 AliGenGeVSim.cxx:770
 AliGenGeVSim.cxx:771
 AliGenGeVSim.cxx:772
 AliGenGeVSim.cxx:773
 AliGenGeVSim.cxx:774
 AliGenGeVSim.cxx:775
 AliGenGeVSim.cxx:776
 AliGenGeVSim.cxx:777
 AliGenGeVSim.cxx:778
 AliGenGeVSim.cxx:779
 AliGenGeVSim.cxx:780
 AliGenGeVSim.cxx:781
 AliGenGeVSim.cxx:782
 AliGenGeVSim.cxx:783
 AliGenGeVSim.cxx:784
 AliGenGeVSim.cxx:785
 AliGenGeVSim.cxx:786
 AliGenGeVSim.cxx:787
 AliGenGeVSim.cxx:788
 AliGenGeVSim.cxx:789
 AliGenGeVSim.cxx:790
 AliGenGeVSim.cxx:791
 AliGenGeVSim.cxx:792
 AliGenGeVSim.cxx:793
 AliGenGeVSim.cxx:794
 AliGenGeVSim.cxx:795
 AliGenGeVSim.cxx:796
 AliGenGeVSim.cxx:797
 AliGenGeVSim.cxx:798
 AliGenGeVSim.cxx:799
 AliGenGeVSim.cxx:800
 AliGenGeVSim.cxx:801
 AliGenGeVSim.cxx:802
 AliGenGeVSim.cxx:803
 AliGenGeVSim.cxx:804
 AliGenGeVSim.cxx:805
 AliGenGeVSim.cxx:806
 AliGenGeVSim.cxx:807
 AliGenGeVSim.cxx:808
 AliGenGeVSim.cxx:809
 AliGenGeVSim.cxx:810
 AliGenGeVSim.cxx:811
 AliGenGeVSim.cxx:812
 AliGenGeVSim.cxx:813
 AliGenGeVSim.cxx:814
 AliGenGeVSim.cxx:815
 AliGenGeVSim.cxx:816
 AliGenGeVSim.cxx:817
 AliGenGeVSim.cxx:818
 AliGenGeVSim.cxx:819
 AliGenGeVSim.cxx:820
 AliGenGeVSim.cxx:821
 AliGenGeVSim.cxx:822
 AliGenGeVSim.cxx:823
 AliGenGeVSim.cxx:824
 AliGenGeVSim.cxx:825
 AliGenGeVSim.cxx:826
 AliGenGeVSim.cxx:827
 AliGenGeVSim.cxx:828
 AliGenGeVSim.cxx:829
 AliGenGeVSim.cxx:830
 AliGenGeVSim.cxx:831
 AliGenGeVSim.cxx:832
 AliGenGeVSim.cxx:833
 AliGenGeVSim.cxx:834
 AliGenGeVSim.cxx:835
 AliGenGeVSim.cxx:836
 AliGenGeVSim.cxx:837
 AliGenGeVSim.cxx:838
 AliGenGeVSim.cxx:839
 AliGenGeVSim.cxx:840
 AliGenGeVSim.cxx:841
 AliGenGeVSim.cxx:842
 AliGenGeVSim.cxx:843
 AliGenGeVSim.cxx:844
 AliGenGeVSim.cxx:845
 AliGenGeVSim.cxx:846
 AliGenGeVSim.cxx:847
 AliGenGeVSim.cxx:848
 AliGenGeVSim.cxx:849
 AliGenGeVSim.cxx:850
 AliGenGeVSim.cxx:851
 AliGenGeVSim.cxx:852
 AliGenGeVSim.cxx:853
 AliGenGeVSim.cxx:854
 AliGenGeVSim.cxx:855
 AliGenGeVSim.cxx:856
 AliGenGeVSim.cxx:857
 AliGenGeVSim.cxx:858
 AliGenGeVSim.cxx:859
 AliGenGeVSim.cxx:860
 AliGenGeVSim.cxx:861
 AliGenGeVSim.cxx:862
 AliGenGeVSim.cxx:863
 AliGenGeVSim.cxx:864
 AliGenGeVSim.cxx:865
 AliGenGeVSim.cxx:866
 AliGenGeVSim.cxx:867
 AliGenGeVSim.cxx:868
 AliGenGeVSim.cxx:869
 AliGenGeVSim.cxx:870
 AliGenGeVSim.cxx:871
 AliGenGeVSim.cxx:872
 AliGenGeVSim.cxx:873
 AliGenGeVSim.cxx:874
 AliGenGeVSim.cxx:875
 AliGenGeVSim.cxx:876
 AliGenGeVSim.cxx:877
 AliGenGeVSim.cxx:878
 AliGenGeVSim.cxx:879
 AliGenGeVSim.cxx:880
 AliGenGeVSim.cxx:881
 AliGenGeVSim.cxx:882
 AliGenGeVSim.cxx:883
 AliGenGeVSim.cxx:884
 AliGenGeVSim.cxx:885
 AliGenGeVSim.cxx:886
 AliGenGeVSim.cxx:887
 AliGenGeVSim.cxx:888
 AliGenGeVSim.cxx:889
 AliGenGeVSim.cxx:890
 AliGenGeVSim.cxx:891
 AliGenGeVSim.cxx:892
 AliGenGeVSim.cxx:893
 AliGenGeVSim.cxx:894
 AliGenGeVSim.cxx:895