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

// Parameterisation of pi and K, eta and pt distributions
// used for the ALICE TDRs.
// eta: according to HIJING (shadowing + quenching)
// pT : according to CDF measurement at 1.8 TeV
// Author: andreas.morsch@cern.ch


//Begin_Html
/*
<img src="picts/AliGeneratorClass.gif">
</pre>
<br clear=left>
<font size=+2 color=red>
<p>The responsible person for this module is
<a href="mailto:andreas.morsch@cern.ch">Andreas Morsch</a>.
</font>
<pre>
*/
//End_Html
//                                                               //
///////////////////////////////////////////////////////////////////

#include <TArrayF.h>
#include <TCanvas.h>
#include <TClonesArray.h>
#include <TDatabasePDG.h>
#include <TF1.h>
#include <TH1.h>
#include <TPDGCode.h>
#include <TParticle.h>
#include <TROOT.h>
#include <TVirtualMC.h>

#include "AliConst.h"
#include "AliDecayer.h"
#include "AliGenEventHeader.h"
#include "AliGenHIJINGpara.h"
#include "AliLog.h"
#include "AliRun.h"

ClassImp(AliGenHIJINGpara)



//_____________________________________________________________________________
static Double_t ptpi(const Double_t *px, const Double_t *)
{
  //
  //     PT-PARAMETERIZATION CDF, PRL 61(88) 1819
  //     POWER LAW FOR PT > 500 MEV
  //     MT SCALING BELOW (T=160 MEV)
  //
  const Double_t kp0    = 1.3;
  const Double_t kxn    = 8.28;
  const Double_t kxlim  = 0.5;
  const Double_t kt     = 0.160;
  const Double_t kxmpi  = 0.139;
  const Double_t kb     = 1.;
  Double_t y, y1, xmpi2, ynorm, a;
  Double_t x = *px;
  //
  y1 = TMath::Power(kp0 / (kp0 + kxlim), kxn);
  xmpi2 = kxmpi * kxmpi;
  ynorm = kb * (TMath::Exp(-sqrt(kxlim * kxlim + xmpi2) / kt ));
  a = ynorm / y1;
  if (x > kxlim)
    y = a * TMath::Power(kp0 / (kp0 + x), kxn);
  else
    y = kb* TMath::Exp(-sqrt(x * x + xmpi2) / kt);
  
  return y*x;
}

//_____________________________________________________________________________
static Double_t ptscal(Double_t pt, Int_t np)
{
    //    SCALING EN MASSE PAR RAPPORT A PTPI
    //     MASS PI,K,ETA,RHO,OMEGA,ETA',PHI
    const Double_t khm[10] = {.13957,.493,.5488,.769,.7826,.958,1.02,0,0,0};
    //     VALUE MESON/PI AT 5 GEV
    const Double_t kfmax[10]={1.,0.3,0.55,1.0,1.0,1.0,1.0,0,0,0};
    np--;
    Double_t f5=TMath::Power(((
	sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
    Double_t fmax2=f5/kfmax[np];
    // PIONS
    Double_t ptpion=100.*ptpi(&pt, (Double_t*) 0);
    Double_t fmtscal=TMath::Power(((
	sqrt(pt*pt+0.018215)+2.)/ (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ 
	fmax2;
    return fmtscal*ptpion;
}

//_____________________________________________________________________________
static Double_t ptka( Double_t *px, Double_t *)
{
    //
    // pt parametrisation for k
    //
    return ptscal(*px,2);
}


//_____________________________________________________________________________
static Double_t etapic( Double_t *py, Double_t *)
{
  //
  // eta parametrisation for pi
  //
    const Double_t ka1    = 4913.;
    const Double_t ka2    = 1819.;
    const Double_t keta1  = 0.22;
    const Double_t keta2  = 3.66;
    const Double_t kdeta1 = 1.47;
    const Double_t kdeta2 = 1.51;
    Double_t y=TMath::Abs(*py);
    //
    Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
    Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
    return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
}

//_____________________________________________________________________________
static Double_t etakac( Double_t *py, Double_t *)
{
    //
    // eta parametrisation for ka
    //
    const Double_t ka1    = 497.6;
    const Double_t ka2    = 215.6;
    const Double_t keta1  = 0.79;
    const Double_t keta2  = 4.09;
    const Double_t kdeta1 = 1.54;
    const Double_t kdeta2 = 1.40;
    Double_t y=TMath::Abs(*py);
    //
    Double_t ex1 = (y-keta1)*(y-keta1)/(2*kdeta1*kdeta1);
    Double_t ex2 = (y-keta2)*(y-keta2)/(2*kdeta2*kdeta2);
    return ka1*TMath::Exp(-ex1)+ka2*TMath::Exp(-ex2);
}

//_____________________________________________________________________________
AliGenHIJINGpara::AliGenHIJINGpara()
    :AliGenerator(),
	fNt(-1),
	fNpartProd(0),
	fPi0Decays(kFALSE),
	fPtWgtPi(0.),
	fPtWgtKa(0.),
	fPtpi(0),
	fPtka(0),
	fETApic(0),
	fETAkac(0),
	fDecayer(0)
{
    //
    // Default constructor
    //
    SetCutVertexZ();
    SetPtRange();
}

//_____________________________________________________________________________
AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
    :AliGenerator(npart),
	fNt(-1),
	fNpartProd(npart),
	fPi0Decays(kFALSE),
	fPtWgtPi(0.),
	fPtWgtKa(0.),
	fPtpi(0),
	fPtka(0),
	fETApic(0),
	fETAkac(0),
	fDecayer(0)
{
  // 
  // Standard constructor
  //
    fName="HIJINGpara";
    fTitle="HIJING Parametrisation Particle Generator";
    SetCutVertexZ();
    SetPtRange();
}

//_____________________________________________________________________________
AliGenHIJINGpara::~AliGenHIJINGpara()
{
  //
  // Standard destructor
  //
    delete fPtpi;
    delete fPtka;
    delete fETApic;
    delete fETAkac;
}

//_____________________________________________________________________________
void AliGenHIJINGpara::Init()
{
  //
  // Initialise the HIJING parametrisation
  //
    Float_t etaMin =-TMath::Log(TMath::Tan(
	TMath::Min((Double_t)fThetaMax/2,TMath::Pi()/2-1.e-10)));
    Float_t etaMax = -TMath::Log(TMath::Tan(
	TMath::Max((Double_t)fThetaMin/2,1.e-10)));
    fPtpi   = new TF1("ptpi",&ptpi,0,20,0);
    gROOT->GetListOfFunctions()->Remove(fPtpi);
    fPtka   = new TF1("ptka",&ptka,0,20,0);
    gROOT->GetListOfFunctions()->Remove(fPtka);
    fPtpi->SetNpx(1000);
    fPtka->SetNpx(1000);
    fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
    gROOT->GetListOfFunctions()->Remove(fETApic);
    fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
    gROOT->GetListOfFunctions()->Remove(fETAkac);

    TF1 etaPic0("etaPic0",&etapic,-7,7,0);
    TF1 etaKac0("etaKac0",&etakac,-7,7,0);

    TF1 ptPic0("ptPic0",&ptpi,0.,15.,0);
    TF1 ptKac0("ptKac0",&ptka,0.,15.,0);

    Float_t intETApi  = etaPic0.Integral(-0.5, 0.5);
    Float_t intETAka  = etaKac0.Integral(-0.5, 0.5);
    Float_t scalePi   = 7316/(intETApi/1.5);
    Float_t scaleKa   =  684/(intETAka/2.0);

//  Fraction of events corresponding to the selected pt-range    
    Float_t intPt    = (0.877*ptPic0.Integral(0, 15)+
			0.123*ptKac0.Integral(0, 15));
    Float_t intPtSel = (0.877*ptPic0.Integral(fPtMin, fPtMax)+
			0.123*ptKac0.Integral(fPtMin, fPtMax));
    Float_t ptFrac   = intPtSel/intPt;

//  Fraction of events corresponding to the selected eta-range    
    Float_t intETASel  = (scalePi*etaPic0.Integral(etaMin, etaMax)+
			  scaleKa*etaKac0.Integral(etaMin, etaMax));
//  Fraction of events corresponding to the selected phi-range    
    Float_t phiFrac    = (fPhiMax-fPhiMin)/2/TMath::Pi();

    
    fParentWeight = (intETASel*ptFrac*phiFrac) / Float_t(fNpart);
    
    if (fAnalog != 0) {
	fPtWgtPi = (fPtMax - fPtMin) / fPtpi->Integral(0., 20.);
	fPtWgtKa = (fPtMax - fPtMin) / fPtka->Integral(0., 20.);
	fParentWeight = (intETASel*phiFrac) / Float_t(fNpart);
    }
    
    
    AliInfo(Form("The number of particles in the selected kinematic region corresponds to %f percent of a full event", 
		 100./ fParentWeight));

// Issue warning message if etaMin or etaMax are outside the alowed range 
// of the parametrization
    if (etaMin < -8.001 || etaMax > 8.001) {
	AliWarning("\nYOU ARE USING THE PARAMETERISATION OUTSIDE ");	
	AliWarning("THE ALLOWED PSEUDORAPIDITY RANGE (-8. - 8.)");	    
	AliWarning(Form("YOUR LIMITS: %f %f \n ", etaMin, etaMax));
    }
//
//
    if (fPi0Decays && TVirtualMC::GetMC())
	fDecayer = TVirtualMC::GetMC()->GetDecayer();

    if (fPi0Decays)
    {
	fDecayer->SetForceDecay(kNeutralPion);
	fDecayer->Init();
    }
    
}


//_____________________________________________________________________________
void AliGenHIJINGpara::Generate()
{
  //
  // Generate one trigger
  //

  
    const Float_t kRaKpic=0.14;
    const Float_t kBorne=1/(1+kRaKpic);
    Float_t polar[3]= {0,0,0};
    //
    const Int_t kPions[3] = {kPi0, kPiPlus, kPiMinus};
    const Int_t kKaons[4] = {kK0Long, kK0Short, kKPlus, kKMinus};
    //
    Float_t origin[3];
    Float_t time;
    Float_t pt, pl, ptot, wgt;
    Float_t phi, theta;
    Float_t p[3];
    Int_t i, part, j;
    //
    TF1 *ptf;
    TF1 *etaf;
    //
    Float_t random[6];
    //
    for (j=0;j<3;j++) origin[j]=fOrigin[j];
    time = fTimeOrigin;

    if(fVertexSmear == kPerEvent) {
	Vertex();
	for (j=0; j < 3; j++) origin[j] = fVertex[j];
	time = fTime;
    } // if kPerEvent
    TArrayF eventVertex;
    eventVertex.Set(3);
    eventVertex[0] = origin[0];
    eventVertex[1] = origin[1];
    eventVertex[2] = origin[2];
    Float_t eventTime = time;
     
    for(i=0;i<fNpart;i++) {
	while(1) {
	    Rndm(random,4);
	    if(random[0]<kBorne) {
		part=kPions[Int_t (random[1]*3)];
		ptf=fPtpi;
		etaf=fETApic;
		wgt = fPtWgtPi;
	    } else {
		part=kKaons[Int_t (random[1]*4)];
		ptf=fPtka;
		etaf=fETAkac;
		wgt = fPtWgtKa;
	    }
	    phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
	    theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
	    if(theta<fThetaMin || theta>fThetaMax) continue;
	 
	    if (fAnalog == 0) { 
		pt = ptf->GetRandom();
	    } else {
		pt = fPtMin + random[3] * (fPtMax - fPtMin);
	    }
	    
	    
	    pl=pt/TMath::Tan(theta);
	    ptot=TMath::Sqrt(pt*pt+pl*pl);
	    if(ptot<fPMin || ptot>fPMax) continue;
	    p[0]=pt*TMath::Cos(phi);
	    p[1]=pt*TMath::Sin(phi);
	    p[2]=pl;
	    if(fVertexSmear==kPerTrack) {
		Rndm(random,6);
		for (j=0;j<3;j++) {
		    origin[j]=fOrigin[j]+fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
			TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
		}
		Rndm(random,2);
		time = fTimeOrigin + fOsigma[2]/TMath::Ccgs()*
		  TMath::Cos(2*random[0]*TMath::Pi())*
		  TMath::Sqrt(-2*TMath::Log(random[1]));
	    }

	    if (fAnalog == 0) { 
		wgt = fParentWeight;
	    } else {
		wgt *= (fParentWeight * ptf->Eval(pt));
	    }
	    
	    
	    if (part == kPi0 && fPi0Decays){
//
//          Decay pi0 if requested
		PushTrack(0,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);
		KeepTrack(fNt);
		DecayPi0(origin, p, time);
	    } else {
      // printf("fNt %d", fNt);
		PushTrack(fTrackIt,-1,part,p,origin,polar,time,kPPrimary,fNt,wgt);

		KeepTrack(fNt);
	    }

	    break;
	}
	SetHighWaterMark(fNt);
    }
//

// Header
    AliGenEventHeader* header = new AliGenEventHeader("HIJINGparam");
// Event Vertex
    header->SetPrimaryVertex(eventVertex);
    header->SetInteractionTime(eventTime);
    header->SetNProduced(fNpartProd);
    if (fContainer) {
      header->SetName(fName);
      fContainer->AddHeader(header);
    } else {
      gAlice->SetGenEventHeader(header);	
    }
}

void AliGenHIJINGpara::SetPtRange(Float_t ptmin, Float_t ptmax) {
    AliGenerator::SetPtRange(ptmin, ptmax);
}

void AliGenHIJINGpara::DecayPi0(Float_t* orig, Float_t * p, Float_t time) 
{
//
//    Decay the pi0
//    and put decay products on the stack
//
    static TClonesArray *particles;
    if(!particles) particles = new TClonesArray("TParticle",1000);
//    
    const Float_t kMass = TDatabasePDG::Instance()->GetParticle(kPi0)->Mass();
    Float_t       e     = TMath::Sqrt(p[0] * p[0] + p[1] * p[1] + p[2] * p[2]+ kMass * kMass);
//
//  Decay the pi0    
    TLorentzVector pmom(p[0], p[1], p[2], e);
    fDecayer->Decay(kPi0, &pmom);
    
//
// Put decay particles on the stack
//
    Float_t polar[3] = {0., 0., 0.};
    Int_t np = fDecayer->ImportParticles(particles);
    fNpartProd += (np-1);
    Int_t nt = 0;    
    for (Int_t i = 1; i < np; i++)
    {
	TParticle* iParticle =  (TParticle *) particles->At(i);
	p[0] = iParticle->Px();
	p[1] = iParticle->Py();
	p[2] = iParticle->Pz();
	Int_t part = iParticle->GetPdgCode();

	PushTrack(fTrackIt, fNt, part, p, orig, polar, time, kPDecay, nt, fParentWeight);
	KeepTrack(nt);
    }
    fNt = nt;
}

void AliGenHIJINGpara::Draw( const char * /*opt*/)
{
    //
    // Draw the pT and y Distributions
    //
     TCanvas *c0 = new TCanvas("c0","Canvas 0",400,10,600,700);
     c0->Divide(2,1);
     c0->cd(1);
     fPtpi->Draw();
     fPtpi->GetHistogram()->SetXTitle("p_{T} (GeV)");     
     c0->cd(2);
     fPtka->Draw();
     fPtka->GetHistogram()->SetXTitle("p_{T} (GeV)");     

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