]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Same class as previously in AliSimpleGen.cxx
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 9 Jun 2000 20:25:01 +0000 (20:25 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 9 Jun 2000 20:25:01 +0000 (20:25 +0000)
All coding rule violations except RS3 corrected (AM)

EVGEN/AliGenBox.cxx [new file with mode: 0644]
EVGEN/AliGenBox.h [new file with mode: 0644]
EVGEN/AliGenFixed.cxx [new file with mode: 0644]
EVGEN/AliGenFixed.h [new file with mode: 0644]
EVGEN/AliGenHIJINGpara.cxx [new file with mode: 0644]
EVGEN/AliGenHIJINGpara.h [new file with mode: 0644]

diff --git a/EVGEN/AliGenBox.cxx b/EVGEN/AliGenBox.cxx
new file mode 100644 (file)
index 0000000..52cc310
--- /dev/null
@@ -0,0 +1,147 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+
+/*
+Old Log:
+Revision 1.8  2000/06/08 13:34:50  fca
+Better control of momentum range in GenBox
+
+Revision 1.7  2000/06/07 16:29:58  fca
+Adding check for pt range in AliGenBox
+
+Revision 1.6  1999/11/03 17:43:20  fca
+New version from G.Martinez & A.Morsch
+
+Revision 1.5  1999/09/29 09:24:14  fca
+Introduction of the Copyright and cvs Log
+*/
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+//    Generate the final state of the interaction as the input   //
+//    to the MonteCarlo                                          //
+//
+//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 "AliGenBox.h"
+#include "AliRun.h"
+#include "AliConst.h"
+#include "AliPDG.h"
+
+ClassImp(AliGenBox)
+
+//_____________________________________________________________________________
+AliGenBox::AliGenBox()
+    :AliGenerator()
+{
+  //
+  // Default constructor
+  //
+  fIpart=0;
+}
+
+//_____________________________________________________________________________
+AliGenBox::AliGenBox(Int_t npart)
+  :AliGenerator(npart)
+{
+  //
+  // Standard constructor
+  //
+  fName="Box";
+  fTitle="Box particle generator";
+  // Generate Proton by default
+  fIpart=kProton;
+}
+
+//_____________________________________________________________________________
+
+void AliGenBox::Generate()
+{
+  //
+  // Generate one trigger
+  //
+  
+    Float_t polar[3]= {0,0,0};
+  //
+    Float_t origin[3];
+    Float_t p[3];
+    Int_t i, j, nt;
+    Double_t pmom, theta, phi, pt;
+    //
+    Float_t random[6];
+  //
+    for (j=0;j<3;j++) origin[j]=fOrigin[j];
+    if(fVertexSmear==perEvent) {
+       gMC->Rndm(random,6);
+       for (j=0;j<3;j++) {
+           origin[j]+=fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
+               TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
+       }
+    }
+    for(i=0;i<fNpart;i++) {
+       gMC->Rndm(random,3);
+       theta=fThetaMin+random[0]*(fThetaMax-fThetaMin);
+       if(TestBit(kMomentumRange)) {
+           pmom=fPMin+random[1]*(fPMax-fPMin);
+           pt=pmom*TMath::Sin(theta);
+       } else {
+
+           pt=fPtMin+random[1]*(fPtMax-fPtMin);
+           pmom=pt/TMath::Sin(theta);
+       }
+       phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
+       p[0] = pt*TMath::Cos(phi);
+       p[1] = pt*TMath::Sin(phi);
+       p[2] = pmom*TMath::Cos(theta);
+
+       if(fVertexSmear==perTrack) {
+           gMC->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]));
+           }
+       }
+       gAlice->SetTrack(fTrackIt,-1,fIpart,p,origin,polar,0,"Primary",nt);
+    }
+}
+
+//_____________________________________________________________________________
+
+void AliGenBox::Init()
+{
+// Initialisation, check consistency of selected ranges
+  if(TestBit(kPtRange)&&TestBit(kMomentumRange)) 
+    Fatal("Init","You should not set the momentum range and the pt range!\n");
+  if((!TestBit(kPtRange))&&(!TestBit(kMomentumRange))) 
+    Fatal("Init","You should set either the momentum or the pt range!\n");
+}
+
diff --git a/EVGEN/AliGenBox.h b/EVGEN/AliGenBox.h
new file mode 100644 (file)
index 0000000..fee9ed6
--- /dev/null
@@ -0,0 +1,26 @@
+#ifndef ALIGENBOX_H
+#define ALIGENBOX_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenerator.h"
+class AliGenBox : public AliGenerator
+{
+ public:
+
+  AliGenBox();
+  AliGenBox(Int_t npart);
+  virtual ~AliGenBox() {}
+  virtual void Generate();
+  virtual void Init();
+  virtual void SetPart(Int_t part) {fIpart=part;}
+protected:
+
+  Int_t fIpart; // Particle type
+
+  ClassDef(AliGenBox,1) // Square box random generator
+};
+
+#endif
diff --git a/EVGEN/AliGenFixed.cxx b/EVGEN/AliGenFixed.cxx
new file mode 100644 (file)
index 0000000..321d195
--- /dev/null
@@ -0,0 +1,95 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+
+///////////////////////////////////////////////////////////////////
+//                                                               //
+//    Generate the final state of the interaction as the input   //
+//    to the MonteCarlo                                          //
+//
+//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 "AliGenFixed.h"
+#include "AliRun.h"
+#include "AliConst.h"
+#include "AliPDG.h"
+#include "TF1.h"
+  
+ClassImp(AliGenFixed)
+
+//_____________________________________________________________________________
+AliGenFixed::AliGenFixed()
+  :AliGenerator()
+{
+  //
+  // Default constructor
+  //
+  fIpart = 0;
+}
+
+//_____________________________________________________________________________
+AliGenFixed::AliGenFixed(Int_t npart)
+  :AliGenerator(npart)
+{
+  //
+  // Standard constructor
+  //
+  fName="Fixed";
+  fTitle="Fixed Particle Generator";
+  // Generate Proton by default
+  fIpart=kProton;
+}
+
+//_____________________________________________________________________________
+void AliGenFixed::Generate()
+{
+  //
+  // Generate one trigger
+  //
+  Float_t polar[3]= {0,0,0};
+  Float_t p[3] = {fPMin*TMath::Cos(fPhiMin)*TMath::Sin(fThetaMin),
+                 fPMin*TMath::Sin(fPhiMin)*TMath::Sin(fThetaMin),
+                 fPMin*TMath::Cos(fThetaMin)};
+  Int_t i, nt;
+  //
+  for(i=0;i<fNpart;i++) {
+    gAlice->SetTrack(fTrackIt,-1,fIpart,p,fOrigin.GetArray(),polar,0,"Primary",nt);
+  }
+}
+  
+//_____________________________________________________________________________
+void AliGenFixed::SetSigma(Float_t sx, Float_t sy, Float_t sz)
+{
+  //
+  // Set the interaction point sigma
+  //
+  printf("Vertex smearing not implemented for fixed generator\n");
+}
diff --git a/EVGEN/AliGenFixed.h b/EVGEN/AliGenFixed.h
new file mode 100644 (file)
index 0000000..ecb43a8
--- /dev/null
@@ -0,0 +1,30 @@
+#ifndef ALIGENFIXED_H
+#define ALIGENFIXED_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenerator.h"
+
+class AliGenFixed : public AliGenerator
+{
+ public:
+  AliGenFixed();
+  AliGenFixed(Int_t npart);
+  virtual ~AliGenFixed() {}
+  virtual void Generate();
+  virtual void Init() {}
+  virtual void SetSigma(Float_t sx, Float_t sy, Float_t sz);
+  virtual void SetMomentum(Float_t pmom) {fPMin=pmom; fPMax=pmom;}
+  virtual void SetPhi(Float_t phi) {fPhiMin=phi*TMath::Pi()/180; fPhiMax=phi*TMath::Pi()/180;}
+  virtual void SetTheta(Float_t theta) {fThetaMin=theta*TMath::Pi()/180; fThetaMax=theta*TMath::Pi()/180;}
+  virtual void SetPart(Int_t part) {fIpart=part;}
+protected:
+
+  Int_t fIpart; // Particle type
+
+  ClassDef(AliGenFixed,1) // Single particle generator
+};
+#endif
diff --git a/EVGEN/AliGenHIJINGpara.cxx b/EVGEN/AliGenHIJINGpara.cxx
new file mode 100644 (file)
index 0000000..6fef731
--- /dev/null
@@ -0,0 +1,299 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+///////////////////////////////////////////////////////////////////
+//                                                               //
+//    Generate the final state of the interaction as the input   //
+//    to the MonteCarlo                                          //
+//
+//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 "AliGenHIJINGpara.h"
+#include "AliRun.h"
+#include "AliConst.h"
+#include "AliPDG.h"
+
+ClassImp(AliGenHIJINGpara)
+
+AliGenHIJINGpara::AliGenHIJINGpara(const AliGenHIJINGpara & para)
+{
+// copy constructor
+}
+
+//_____________________________________________________________________________
+static Double_t ptpi(Double_t *px, 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()
+{
+    //
+    // Default constructor
+    //
+    fPtpi = 0;
+    fPtka = 0;
+    fETApic = 0;
+    fETAkac = 0;
+}
+
+//_____________________________________________________________________________
+AliGenHIJINGpara::AliGenHIJINGpara(Int_t npart)
+  :AliGenerator(npart)
+{
+  // 
+  // Standard constructor
+  //
+    fName="HIGINGpara";
+    fTitle="HIJING Parametrisation Particle Generator";
+    fPtpi = 0;
+    fPtka = 0;
+    fETApic = 0;
+    fETAkac = 0;
+}
+
+//_____________________________________________________________________________
+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);
+    fPtka = new TF1("ptka",&ptka,0,20,0);
+    fETApic = new TF1("etapic",&etapic,etaMin,etaMax,0);
+    fETAkac = new TF1("etakac",&etakac,etaMin,etaMax,0);
+    TF1 *etaPic0 = new TF1("etapic",&etapic,-7,7,0);
+    TF1 *etaKac0 = new TF1("etakac",&etakac,-7,7,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);
+    
+    Float_t intPt  = (0.877*etaPic0->Integral(0, 15)+
+                     0.123*etaKac0->Integral(0, 15));
+    Float_t intPtSel = (0.877*etaPic0->Integral(fPtMin, fPtMax)+
+                       0.123*etaKac0->Integral(fPtMin, fPtMax));
+    Float_t ptFrac = intPtSel/intPt;
+    
+    
+    Float_t intETASel  = (scalePi*etaPic0->Integral(etaMin, etaMax)+
+                         scaleKa*etaKac0->Integral(etaMin, etaMax));
+    Float_t phiFrac = (fPhiMax-fPhiMin)/2/TMath::Pi();
+    fParentWeight = Float_t(fNpart)/intETASel*ptFrac*phiFrac;
+    
+    printf("\n The number of particles in the selected kinematic region corresponds to %f percent of a full event\n ", 100.*fParentWeight);
+    
+}
+
+//_____________________________________________________________________________
+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 pt, pl, ptot;
+    Float_t phi, theta;
+    Float_t p[3];
+    Int_t i, part, nt, j;
+    //
+    TF1 *ptf;
+    TF1 *etaf;
+    //
+    Float_t random[6];
+    //
+    for (j=0;j<3;j++) origin[j]=fOrigin[j];
+    if(fVertexSmear==perEvent) {
+       gMC->Rndm(random,6);
+       for (j=0;j<3;j++) {
+           origin[j]+=fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
+               TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
+       }
+    }
+    for(i=0;i<fNpart;i++) {
+       while(1) {
+           gMC->Rndm(random,3);
+           if(random[0]<kBorne) {
+               part=kPions[Int_t (random[1]*3)];
+               ptf=fPtpi;
+             etaf=fETApic;
+           } else {
+               part=kKaons[Int_t (random[1]*4)];
+               ptf=fPtka;
+               etaf=fETAkac;
+           }
+           phi=fPhiMin+random[2]*(fPhiMax-fPhiMin);
+           theta=2*TMath::ATan(TMath::Exp(-etaf->GetRandom()));
+           if(theta<fThetaMin || theta>fThetaMax) continue;
+           pt=ptf->GetRandom();
+           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==perTrack) {
+               gMC->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]));
+               }
+           }
+           gAlice->SetTrack(fTrackIt,-1,part,p,origin,polar,0,"Primary",nt,fParentWeight);
+           break;
+       }
+    }
+}
+
+AliGenHIJINGpara& AliGenHIJINGpara::operator=(const  AliGenHIJINGpara& rhs)
+{
+// Assignment operator
+    return *this;
+}
+
+
diff --git a/EVGEN/AliGenHIJINGpara.h b/EVGEN/AliGenHIJINGpara.h
new file mode 100644 (file)
index 0000000..effa0fa
--- /dev/null
@@ -0,0 +1,42 @@
+#ifndef ALIGENHIJINGPARA_H
+#define ALIGENHIJINGPARA_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenerator.h"
+#include "TF1.h"
+
+class AliGenHIJINGpara : public AliGenerator
+{
+ public:
+
+  AliGenHIJINGpara();
+  AliGenHIJINGpara(Int_t npart);
+  AliGenHIJINGpara(const AliGenHIJINGpara &HIJINGpara);
+     
+  virtual ~AliGenHIJINGpara();
+  virtual void Generate();
+  virtual void Init();
+  AliGenHIJINGpara & operator=(const AliGenHIJINGpara & rhs);
+ protected:
+
+  TF1* fPtpi; // Parametrised pt distribution for pi
+  TF1* fPtka; // Parametrised pt distribution for ka
+  TF1* fETApic; // Parametrised eta distribution for pi
+  TF1* fETAkac; // Parametrised eta distribution fro ka
+
+  ClassDef(AliGenHIJINGpara,1) // Hijing parametrisation generator
+};
+#endif
+
+
+
+
+
+
+
+
+
+