1 #ifndef __StarGenerator_h__
2 #define __StarGenerator_h__
6 #include "TDatabasePDG.h"
7 #include "StarCallf77.h"
8 #include "TLorentzVector.h"
10 #include "StarGenerator/EVENT/StarGenEvent.h"
11 #include "StarGenerator/EVENT/StarGenStats.h"
12 #include "StarGenerator/UTIL/StarParticleData.h"
51 virtual Int_t
Init() = 0;
80 Int_t PreGenerateHook();
82 Int_t PostGenerateHook();
89 void Clear(
const Option_t *opts=
"" );
104 void SetBlue(
const Char_t *b );
106 void SetYell(
const Char_t *y );
110 void SetFrame(
const Char_t *frame,
const Double_t val );
115 void SetFrame(
const Char_t *frame,
const Double_t *pblue,
const Double_t *pyell );
147 void SetInputTree ( TTree *tree,
const Char_t *bname=0 );
149 void SetInputFile (
const Char_t *name,
const Char_t *treename=
"genevents",
const Char_t *bname=0 );
163 virtual const char *GetCVS()
const
164 {
static const char cvs[]=
"Tag $Name: $ $Id: StarGenerator.h,v 1.6 2015/04/24 18:21:11 jwebb Exp $ built " __DATE__
" " __TIME__ ;
return cvs;}
172 void SetId( Int_t
id ) { mId = id; }
176 Int_t mNumberOfParticles;
Double_t mRootS
CMS energy or incident beam momentum for fixed target collisions.
TString mLibrary
Concrete event generator library.
virtual Int_t Generate()=0
void SetOutputTree(TTree *tree)
TString mYell
Name of the yellow beam particle (-z)
Bool_t IsPileup()
Returns true if this is a pileup generator.
Int_t IOmode()
Returns the ID of this event generator.
StarGenEvent * Event()
Retrieves the event record.
void SetFrame(const Char_t *frame, const Double_t val)
Double_t mVertex[3]
Position of the interaction vertex.
Int_t Make()
If called, make simply aliases to the Generate() function.
void SetImpact(Double_t bmin, Double_t bmax)
Set the minimum and maximum impact parameters for the collision (HI generators)
Bool_t mIsPileup
Flag which indicates whether this is a pileup simulation.
TLorentzVector mYellMomentum
4-Momentum of the yellow beam
Int_t mIOMode
IO flag. 0=no IO, 1=write, 2=read.
virtual void FillEP(StarGenEvent *event)
(Optional) Method to fill a DIS event
ABC for defining event generator interfaces.
void Clear(const Option_t *opts="")
Clear the event.
void SetInputTree(TTree *tree, const Char_t *bname=0)
StarGenEvent * mEvent
Generated event.
Double_t mImpactMin
Minimum impact parameter in a HI collision.
virtual void FillEA(StarGenEvent *event)
(Optional) Method to fill a eA event
TTree * mOutputTree
Pointer to the tree for writing out events.
Interface to PDG information.
Int_t InitRun(Int_t runnumber)
Developer may not provide run initialization.
Double_t mImpactMax
Maximum impact parameter in a HI collision.
void SetBlue(const Char_t *b)
Sets the particle species for the blue beam.
Double_t mDirect
Direction (+1 = W, -1 = E) of the beam in fixted target mode.
Int_t GetNumberOfParticles()
Returns the number of particles in the event.
TLorentzVector mBlueMomentum
4-Momentum of the blue beam
virtual Int_t PostGenerate()
Developers may provide a post-generate method which will execute after Generate().
Double_t mSpread[3]
Width of the interaction vertex.
Double_t mPileup
Indicates the probability that this generator should be run in pileup mode.
const Double_t cm
If an event generator measures distances in cm, then multiply by cm.
End of run statistics for event generators.
TTree * mInputTree
Pointer to the tree for reading in events.
const Double_t mm
If an event generator measures distances in mm, then multiply by mm.
Base class for event records.
virtual void FillUser(StarGenEvent *event)
(Optional) Method to fill a user-defined event
void SetPileup(Bool_t p, Double_t prob=0.10)
virtual StarGenStats Stats()
Create and retrieve end-of-run statistics.
Main steering class for event generation.
Double_t mBlueMass
Mass of the blue beam.
void SetYell(const Char_t *y)
Sets the particle species for the yellow beam.
Double_t mYellMass
Mass of the yellow beam.
StarParticleData * mParticleDb
Database for particle information (based on TDatabasePDG)
TString mFrame
Frame of the collision, i.e. CMS, FIXT, 3MOM, 4MOM, 5MOM.
TString mBlue
Name of the blue beam particle (+z)
virtual void FillPP(StarGenEvent *event)
(Optional) Method to fill a PP event
virtual Int_t PreGenerate()
Developers may provide a pre-generate method which will execute before Generate().
virtual void FillAA(StarGenEvent *event)
(Optional) Method to fill a AA event
Double_t GetPileup()
Returns the pileup probability.