Geant4 User's Guide
For Application Developers Toolkit Fundamentals |
3.6 Event Generator Interface
3.6.1 Structure of a primary event
Primary vertex and primary particle
The G4Event class object should have a set of primary particles when it is sent to G4EventManager via processOneEvent() method. It is the mandate of your G4VUserPrimaryGeneratorAction concrete class to send primary particles to the G4Event object.
The G4PrimaryParticle class represents a primary particle with which Geant4 starts simulating an event. This class object has information on particle type and its three momenta. The positional and time information of primary particle(s) are stored in the G4PrimaryVertex class object and, thus, this class object can have one or more G4PrimaryParticle class objects which share the same vertex. As shown in Fig.?.?, primary vertexes and primary particles are associated with the G4Event object by a form of linked list.
A concrete class of G4VPrimaryGenerator, the G4PrimaryParticle object is constructed with either a pointer to G4ParticleDefinition or an integer number which represents P.D.G. particle code. For the case of some artificial particles, e.g., geantino, optical photon, etc., or exotic nuclear fragments, which the P.D.G. particle code does not cover, the G4PrimaryParticle should be constructed by G4ParticleDefinition pointer. On the other hand, elementary particles with very short life time, e.g., weak bosons, or quarks/gluons, can be instantiated as G4PrimaryParticle objects using the P.D.G. particle code. It should be noted that, even though primary particles with such a very short life time are defined, Geant4 will simulate only the particles which are defined as G4ParticleDefinition class objects. Other primary particles will be simply ignored by G4EventManager. But it may still be useful to construct such "intermediate" particles for recording the origin of the primary event.
Forced decay channel
The G4PrimaryParticle class object can have a list of its daughter particles. If the parent particle is an "intermediate" particle, which Geant4 does not have a corresponding G4ParticleDefinition, this parent particle is ignored and daughters are assumed to start from the vertex with which their parent is associated. For example, a Z boson is associated with a vertex and it has positive and negative muons as its daughters, these muons will start from that vertex.
There are some kinds of particles which should fly some reasonable distances and, thus, should be simulated by Geant4, but you still want to follow the decay channel generated by an event generator. A typical case of these particles is B meson. Even for the case of a primary particle which has a corresponding G4ParticleDefinition, it can have daughter primary particles. Geant4 will trace the parent particle until it comes to decay, obeying multiple scattering, ionization loss, rotation with the magnetic field, etc. according to its particle type. When the parent comes to decay, instead of randomly choosing its decay channel, it follows the "pre-assigned" decay channel. To conserve the energy and the momentum of the parent, daughters will be Lorentz transformed according to their parent's frame.
Unfortunately, almost all event generators presently in use, commonly are written in FORTRAN. For Geant4, it was decided to not link with any FORTRAN program or library, even though the C++ language syntax itself allows such a link. Linking to a FORTRAN package might be convenient in some cases, but we will lose many advantages of object-oriented features of C++, such as robustness. Instead, Geant4 provides an ASCII file interface for such event generators.
G4HEPEvtInterface is one of G4VPrimaryGenerator concrete class and thus it can be used in your G4VUserPrimaryGeneratorAction concrete class. G4HEPEvtInterface reads an ASCII file produced by an event generator and reproduces G4PrimaryParticle objects associated with a G4PrimaryVertex object. It reproduces a full production chain of the event generator, starting with primary quarks, etc. In other words, G4HEPEvtInterface converts information stored in the /HEPEVT/ common block to an object-oriented data structure. Because the /HEPEVT/ common block is commonly used by almost all event generators written in FORTRAN, G4HEPEvtInterface can interface to almost all event generators currently used in the HEP community. The constructor of G4HEPEvtInterface takes the file name. Source listing 3.6.1 shows an example how to use G4HEPEvtInterface. Note that an event generator is not assumed to give a place of the primary particles, the interaction point must be set before invoking GeneratePrimaryVertex() method.
#ifndef ExN04PrimaryGeneratorAction_h #define ExN04PrimaryGeneratorAction_h 1 #include "G4VUserPrimaryGeneratorAction.hh" #include "globals.hh" class G4VPrimaryGenerator; class G4Event; class ExN04PrimaryGeneratorAction : public G4VUserPrimaryGeneratorAction { public: ExN04PrimaryGeneratorAction(); ~ExN04PrimaryGeneratorAction(); public: void GeneratePrimaries(G4Event* anEvent); private: G4VPrimaryGenerator* HEPEvt; }; #endif #include "ExN04PrimaryGeneratorAction.hh" #include "G4Event.hh" #include "G4HEPEvtInterface.hh" ExN04PrimaryGeneratorAction::ExN04PrimaryGeneratorAction() { HEPEvt = new G4HEPEvtInterface("pythia_event.data"); } ExN04PrimaryGeneratorAction::~ExN04PrimaryGeneratorAction() { delete HEPEvt; } void ExN04PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent) { HEPEvt->SetParticlePosition(G4ThreeVector(0.*cm,0.*cm,0.*cm)); HEPEvt->GeneratePrimaryVertex(anEvent); } |
Source listing 3.6.1 An example code for G4HEPEvtInterface. |
Format of the ASCII file
An ASCII file, which will be fed by G4HEPEvtInterface should have the following format.
*********************************************************** SUBROUTINE HEP2G4 * * Convert /HEPEVT/ event structure to an ASCII file * to be fed by G4HEPEvtInterface * *********************************************************** PARAMETER (NMXHEP=2000) COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(NMXHEP),IDHEP(NMXHEP), >JMOHEP(2,NMXHEP),JDAHEP(2,NMXHEP),PHEP(5,NMXHEP),VHEP(4,NMXHEP) DOUBLE PRECISION PHEP,VHEP * WRITE(6,*) NHEP DO IHEP=1,NHEP WRITE(6,10) > ISTHEP(IHEP),IDHEP(IHEP),JDAHEP(1,IHEP),JDAHEP(2,IHEP), > PHEP(1,IHEP),PHEP(2,IHEP),PHEP(3,IHEP),PHEP(5,IHEP) 10 FORMAT(4I10,4(1X,D15.8)) ENDDO * RETURN END |
Source listing 3.6.2 A FORTRAN example using the /HEPEVT/ common. |
Future interface to the new generation generators
Several activities have already been started for developing object-oriented event generators. Such new generators can be easily linked and used with a Geant4 based simulation. Furthermore, we need not distinguish a primary generator from the physics processes used in Geant4. Future generators can be a kind of physics process plugged-in by inheriting G4VProcess.
One possible use is the following. Within an event,
a G4HEPEvtInterface class object instantiated with a
minimum bias event file is accessed 20 times and another
G4HEPEvtInterface class object instantiated with a
signal event file is accessed once. Thus, this event represents
a typical signal event of LHC overlapping 20 minimum bias
events. It should be noted that a simulation of event
overlapping can be done by merging hits and/or digits
associated with several events, and these events can be
simulated independently. Digitization over multiple events
will be mentioned in Section 4.5.