| 1 | <html><head><meta http-equiv="Content-Type" content="text/html; charset=ISO-8859-1"><title>3.6. Event Generator Interface</title><link rel="stylesheet" href="../xml/XSLCustomizationLayer/G4HTMLStylesheet.css" type="text/css"><meta name="generator" content="DocBook XSL Stylesheets V1.71.1"><link rel="start" href="index.html" title="Geant4 User's Guide for Application Developers"><link rel="up" href="ch03.html" title="Chapter 3. Toolkit Fundamentals"><link rel="prev" href="ch03s05.html" title="3.5. Event"><link rel="next" href="ch03s07.html" title="3.7. Event Biasing Techniques"><script language="JavaScript">
|
|---|
| 2 | function remote_win(fName)
|
|---|
| 3 | {
|
|---|
| 4 | var url = "AllResources/Detector/geometry.src/" + fName;
|
|---|
| 5 | RemoteWin=window.open(url,"","resizable=no,toolbar=0,location=0,directories=0,status=0,menubar=0,scrollbars=0,copyhistory=0,width=520,height=520")
|
|---|
| 6 | RemoteWin.creator=self
|
|---|
| 7 | }
|
|---|
| 8 | </script></head><body bgcolor="white" text="black" link="#0000FF" vlink="#840084" alink="#0000FF"><div class="navheader"><table width="100%" summary="Navigation header"><tr><th colspan="3" align="center">3.6.
|
|---|
| 9 | Event Generator Interface
|
|---|
| 10 | </th></tr><tr><td width="20%" align="left"><a accesskey="p" href="ch03s05.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><th width="60%" align="center">Chapter 3.
|
|---|
| 11 | Toolkit Fundamentals
|
|---|
| 12 | </th><td width="20%" align="right"> <a accesskey="n" href="ch03s07.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr></table><hr></div><div class="sect1" lang="en"><div class="titlepage"><div><div><h2 class="title" style="clear: both"><a name="sect.EventGen"></a>3.6.
|
|---|
| 13 | Event Generator Interface
|
|---|
| 14 | </h2></div></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EventGen.StructPriEvt"></a>3.6.1.
|
|---|
| 15 | Structure of a primary event
|
|---|
| 16 | </h3></div></div></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EventGen.StructPriEvt.VtxPart"></a>3.6.1.1.
|
|---|
| 17 | Primary vertex and primary particle
|
|---|
| 18 | </h4></div></div></div><p>
|
|---|
| 19 | The <span class="emphasis"><em>G4Event</em></span> class object should have a set of primary
|
|---|
| 20 | particles when it is sent to <span class="emphasis"><em>G4EventManager</em></span> via
|
|---|
| 21 | <code class="literal">processOneEvent()</code> method. It is the mandate of your
|
|---|
| 22 | <span class="emphasis"><em>G4VUserPrimaryGeneratorAction</em></span> concrete class to send primary
|
|---|
| 23 | particles to the <span class="emphasis"><em>G4Event</em></span> object.
|
|---|
| 24 | </p><p>
|
|---|
| 25 | The <span class="emphasis"><em>G4PrimaryParticle</em></span> class represents a primary particle
|
|---|
| 26 | with which Geant4 starts simulating an event. This class object has
|
|---|
| 27 | information on particle type and its three momenta. The positional
|
|---|
| 28 | and time information of primary particle(s) are stored in the
|
|---|
| 29 | <span class="emphasis"><em>G4PrimaryVertex</em></span> class object and, thus, this class object
|
|---|
| 30 | can have one or more <span class="emphasis"><em>G4PrimaryParticle</em></span> class objects which
|
|---|
| 31 | share the same vertex. As shown in Fig.?.?, primary vertexes and
|
|---|
| 32 | primary particles are associated with the <span class="emphasis"><em>G4Event</em></span> object by
|
|---|
| 33 | a form of linked list.
|
|---|
| 34 | </p><p>
|
|---|
| 35 | A concrete class of <span class="emphasis"><em>G4VPrimaryGenerator</em></span>, the
|
|---|
| 36 | <span class="emphasis"><em>G4PrimaryParticle</em></span> object is constructed with either a
|
|---|
| 37 | pointer to <span class="emphasis"><em>G4ParticleDefinition</em></span> or an integer number which
|
|---|
| 38 | represents P.D.G. particle code. For the case of some artificial
|
|---|
| 39 | particles, e.g., geantino, optical photon, etc., or exotic nuclear
|
|---|
| 40 | fragments, which the P.D.G. particle code does not cover, the
|
|---|
| 41 | <span class="emphasis"><em>G4PrimaryParticle</em></span> should be constructed by
|
|---|
| 42 | <span class="emphasis"><em>G4ParticleDefinition</em></span> pointer. On the other hand, elementary
|
|---|
| 43 | particles with very short life time, e.g., weak bosons, or
|
|---|
| 44 | quarks/gluons, can be instantiated as <span class="emphasis"><em>G4PrimaryParticle</em></span>
|
|---|
| 45 | objects using the P.D.G. particle code. It should be noted that,
|
|---|
| 46 | even though primary particles with such a very short life time are
|
|---|
| 47 | defined, Geant4 will simulate only the particles which are defined
|
|---|
| 48 | as <span class="emphasis"><em>G4ParticleDefinition</em></span> class objects. Other primary
|
|---|
| 49 | particles will be simply ignored by <span class="emphasis"><em>G4EventManager</em></span>. But it
|
|---|
| 50 | may still be useful to construct such "intermediate" particles for
|
|---|
| 51 | recording the origin of the primary event.
|
|---|
| 52 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EventGen.StructPriEvt.ForcedDecay"></a>3.6.1.2.
|
|---|
| 53 | Forced decay channel
|
|---|
| 54 | </h4></div></div></div><p>
|
|---|
| 55 | The <span class="emphasis"><em>G4PrimaryParticle</em></span> class object can have a list of its
|
|---|
| 56 | daughter particles. If the parent particle is an "intermediate"
|
|---|
| 57 | particle, which Geant4 does not have a corresponding
|
|---|
| 58 | <span class="emphasis"><em>G4ParticleDefinition</em></span>, this parent particle is ignored and
|
|---|
| 59 | daughters are assumed to start from the vertex with which their
|
|---|
| 60 | parent is associated. For example, a Z boson is associated with a
|
|---|
| 61 | vertex and it has positive and negative muons as its daughters,
|
|---|
| 62 | these muons will start from that vertex.
|
|---|
| 63 | </p><p>
|
|---|
| 64 | There are some kinds of particles which should fly some
|
|---|
| 65 | reasonable distances and, thus, should be simulated by Geant4, but
|
|---|
| 66 | you still want to follow the decay channel generated by an event
|
|---|
| 67 | generator. A typical case of these particles is B meson. Even for
|
|---|
| 68 | the case of a primary particle which has a corresponding
|
|---|
| 69 | <span class="emphasis"><em>G4ParticleDefinition</em></span>, it can have daughter primary
|
|---|
| 70 | particles. Geant4 will trace the parent particle until it comes to
|
|---|
| 71 | decay, obeying multiple scattering, ionization loss, rotation with
|
|---|
| 72 | the magnetic field, etc. according to its particle type. When the
|
|---|
| 73 | parent comes to decay, instead of randomly choosing its decay
|
|---|
| 74 | channel, it follows the "pre-assigned" decay channel. To conserve
|
|---|
| 75 | the energy and the momentum of the parent, daughters will be
|
|---|
| 76 | Lorentz transformed according to their parent's frame.
|
|---|
| 77 | </p></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EventGen.InterPriGen"></a>3.6.2.
|
|---|
| 78 | Interface to a primary generator
|
|---|
| 79 | </h3></div></div></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EventGen.InterPriGen.HEPEvt"></a>3.6.2.1.
|
|---|
| 80 | G4HEPEvtInterface
|
|---|
| 81 | </h4></div></div></div><p>
|
|---|
| 82 | Unfortunately, almost all event generators presently in use,
|
|---|
| 83 | commonly are written in FORTRAN. For Geant4, it was decided to not
|
|---|
| 84 | link with any FORTRAN program or library, even though the C++
|
|---|
| 85 | language syntax itself allows such a link. Linking to a FORTRAN
|
|---|
| 86 | package might be convenient in some cases, but we will lose many
|
|---|
| 87 | advantages of object-oriented features of C++, such as robustness.
|
|---|
| 88 | Instead, Geant4 provides an ASCII file interface for such event
|
|---|
| 89 | generators.
|
|---|
| 90 | </p><p>
|
|---|
| 91 | <span class="emphasis"><em>G4HEPEvtInterface</em></span> is one of <span class="emphasis"><em>G4VPrimaryGenerator</em></span>
|
|---|
| 92 | concrete class and thus it can be used in your
|
|---|
| 93 | <span class="emphasis"><em>G4VUserPrimaryGeneratorAction</em></span> concrete class.
|
|---|
| 94 | <span class="emphasis"><em>G4HEPEvtInterface</em></span> reads an ASCII file produced by an event
|
|---|
| 95 | generator and reproduces <span class="emphasis"><em>G4PrimaryParticle</em></span> objects
|
|---|
| 96 | associated with a <span class="emphasis"><em>G4PrimaryVertex</em></span> object. It reproduces a
|
|---|
| 97 | full production chain of the event generator, starting with primary
|
|---|
| 98 | quarks, etc. In other words, <span class="emphasis"><em>G4HEPEvtInterface</em></span> converts
|
|---|
| 99 | information stored in the <code class="literal">/HEPEVT/</code> common block to an
|
|---|
| 100 | object-oriented data structure. Because the <code class="literal">/HEPEVT/</code>
|
|---|
| 101 | common block is commonly used by almost all event generators
|
|---|
| 102 | written in FORTRAN, <span class="emphasis"><em>G4HEPEvtInterface</em></span> can interface to
|
|---|
| 103 | almost all event generators currently used in the HEP community.
|
|---|
| 104 | The constructor of <span class="emphasis"><em>G4HEPEvtInterface</em></span> takes the file name.
|
|---|
| 105 | <a href="ch03s06.html#programlist_EventGen_1" title="Example 3.3.
|
|---|
| 106 | An example code for G4HEPEvtInterface
|
|---|
| 107 | ">Example 3.3</a> shows an example how to use
|
|---|
| 108 | <span class="emphasis"><em>G4HEPEvtInterface</em></span>. Note that an event generator is not
|
|---|
| 109 | assumed to give a place of the primary particles, the interaction
|
|---|
| 110 | point must be set before invoking <span class="emphasis"><em>GeneratePrimaryVertex()</em></span>
|
|---|
| 111 | method.
|
|---|
| 112 | </p><div class="example"><a name="programlist_EventGen_1"></a><p class="title"><b>Example 3.3.
|
|---|
| 113 | An example code for <span class="emphasis"><em>G4HEPEvtInterface</em></span>
|
|---|
| 114 | </b></p><div class="example-contents"><pre class="programlisting">
|
|---|
| 115 | #ifndef ExN04PrimaryGeneratorAction_h
|
|---|
| 116 | #define ExN04PrimaryGeneratorAction_h 1
|
|---|
| 117 |
|
|---|
| 118 | #include "G4VUserPrimaryGeneratorAction.hh"
|
|---|
| 119 | #include "globals.hh"
|
|---|
| 120 |
|
|---|
| 121 | class G4VPrimaryGenerator;
|
|---|
| 122 | class G4Event;
|
|---|
| 123 |
|
|---|
| 124 | class ExN04PrimaryGeneratorAction : public G4VUserPrimaryGeneratorAction
|
|---|
| 125 | {
|
|---|
| 126 | public:
|
|---|
| 127 | ExN04PrimaryGeneratorAction();
|
|---|
| 128 | ~ExN04PrimaryGeneratorAction();
|
|---|
| 129 |
|
|---|
| 130 | public:
|
|---|
| 131 | void GeneratePrimaries(G4Event* anEvent);
|
|---|
| 132 |
|
|---|
| 133 | private:
|
|---|
| 134 | G4VPrimaryGenerator* HEPEvt;
|
|---|
| 135 | };
|
|---|
| 136 |
|
|---|
| 137 | #endif
|
|---|
| 138 |
|
|---|
| 139 |
|
|---|
| 140 | #include "ExN04PrimaryGeneratorAction.hh"
|
|---|
| 141 |
|
|---|
| 142 | #include "G4Event.hh"
|
|---|
| 143 | #include "G4HEPEvtInterface.hh"
|
|---|
| 144 |
|
|---|
| 145 | ExN04PrimaryGeneratorAction::ExN04PrimaryGeneratorAction()
|
|---|
| 146 | {
|
|---|
| 147 | HEPEvt = new G4HEPEvtInterface("pythia_event.data");
|
|---|
| 148 | }
|
|---|
| 149 |
|
|---|
| 150 | ExN04PrimaryGeneratorAction::~ExN04PrimaryGeneratorAction()
|
|---|
| 151 | {
|
|---|
| 152 | delete HEPEvt;
|
|---|
| 153 | }
|
|---|
| 154 |
|
|---|
| 155 | void ExN04PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
|
|---|
| 156 | {
|
|---|
| 157 | HEPEvt->SetParticlePosition(G4ThreeVector(0.*cm,0.*cm,0.*cm));
|
|---|
| 158 | HEPEvt->GeneratePrimaryVertex(anEvent);
|
|---|
| 159 | }
|
|---|
| 160 | </pre></div></div><br class="example-break"></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EventGen.InterPriGen.FormASCII"></a>3.6.2.2.
|
|---|
| 161 | Format of the ASCII file
|
|---|
| 162 | </h4></div></div></div><p>
|
|---|
| 163 | An ASCII file, which will be fed by <span class="emphasis"><em>G4HEPEvtInterface</em></span>
|
|---|
| 164 | should have the following format.
|
|---|
| 165 |
|
|---|
| 166 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 167 | The first line of each primary event should be an integer which
|
|---|
| 168 | represents the number of the following lines of primary
|
|---|
| 169 | particles.
|
|---|
| 170 | </p></li><li><p>
|
|---|
| 171 | Each line in an event corresponds to a particle in the
|
|---|
| 172 | <code class="literal">/HEPEVT/</code> common. Each line has <code class="literal">ISTHEP, IDHEP,
|
|---|
| 173 | JDAHEP(1), JDAHEP(2), PHEP(1), PHEP(2), PHEP(3), PHEP(5).</code>
|
|---|
| 174 | Refer to the <code class="literal">/HEPEVT/</code> manual for the meanings of these
|
|---|
| 175 | variables.
|
|---|
| 176 | </p></li></ul></div><p>
|
|---|
| 177 | </p><p>
|
|---|
| 178 | <a href="ch03s06.html#programlist_EventGen_2" title="Example 3.4.
|
|---|
| 179 | A FORTRAN example using the /HEPEVT/ common.
|
|---|
| 180 | ">Example 3.4</a> shows an example FORTRAN code to
|
|---|
| 181 | generate an ASCII file.
|
|---|
| 182 | </p><div class="example"><a name="programlist_EventGen_2"></a><p class="title"><b>Example 3.4.
|
|---|
| 183 | A FORTRAN example using the <code class="literal">/HEPEVT/</code> common.
|
|---|
| 184 | </b></p><div class="example-contents"><pre class="programlisting">
|
|---|
| 185 | ***********************************************************
|
|---|
| 186 | SUBROUTINE HEP2G4
|
|---|
| 187 | *
|
|---|
| 188 | * Convert /HEPEVT/ event structure to an ASCII file
|
|---|
| 189 | * to be fed by G4HEPEvtInterface
|
|---|
| 190 | *
|
|---|
| 191 | ***********************************************************
|
|---|
| 192 | PARAMETER (NMXHEP=2000)
|
|---|
| 193 | COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(NMXHEP),IDHEP(NMXHEP),
|
|---|
| 194 | >JMOHEP(2,NMXHEP),JDAHEP(2,NMXHEP),PHEP(5,NMXHEP),VHEP(4,NMXHEP)
|
|---|
| 195 | DOUBLE PRECISION PHEP,VHEP
|
|---|
| 196 | *
|
|---|
| 197 | WRITE(6,*) NHEP
|
|---|
| 198 | DO IHEP=1,NHEP
|
|---|
| 199 | WRITE(6,10)
|
|---|
| 200 | > ISTHEP(IHEP),IDHEP(IHEP),JDAHEP(1,IHEP),JDAHEP(2,IHEP),
|
|---|
| 201 | > PHEP(1,IHEP),PHEP(2,IHEP),PHEP(3,IHEP),PHEP(5,IHEP)
|
|---|
| 202 | 10 FORMAT(4I10,4(1X,D15.8))
|
|---|
| 203 | ENDDO
|
|---|
| 204 | *
|
|---|
| 205 | RETURN
|
|---|
| 206 | END
|
|---|
| 207 | </pre></div></div><br class="example-break"></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EventGen.InterPriGen.FutureInt"></a>3.6.2.3.
|
|---|
| 208 | Future interface to the new generation generators
|
|---|
| 209 | </h4></div></div></div><p>
|
|---|
| 210 | Several activities have already been started for developing
|
|---|
| 211 | object-oriented event generators. Such new generators can be easily
|
|---|
| 212 | linked and used with a Geant4 based simulation. Furthermore, we
|
|---|
| 213 | need not distinguish a primary generator from the physics processes
|
|---|
| 214 | used in Geant4. Future generators can be a kind of physics process
|
|---|
| 215 | plugged-in by inheriting <span class="emphasis"><em>G4VProcess</em></span>.
|
|---|
| 216 | </p></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EventGen.EvtOverlap"></a>3.6.3.
|
|---|
| 217 | Event overlap using multiple generators
|
|---|
| 218 | </h3></div></div></div><p>
|
|---|
| 219 | Your <span class="emphasis"><em>G4VUserPrimaryGeneratorAction</em></span> concrete class can have
|
|---|
| 220 | more than one <span class="emphasis"><em>G4VPrimaryGenerator</em></span> concrete class. Each
|
|---|
| 221 | <span class="emphasis"><em>G4VPrimaryGenerator</em></span> concrete class can be accessed more than
|
|---|
| 222 | once per event. Using these class objects, one event can have more
|
|---|
| 223 | than one primary event.
|
|---|
| 224 | </p><p>
|
|---|
| 225 | One possible use is the following. Within an event, a
|
|---|
| 226 | <span class="emphasis"><em>G4HEPEvtInterface</em></span> class object instantiated with a minimum
|
|---|
| 227 | bias event file is accessed 20 times and another
|
|---|
| 228 | <span class="emphasis"><em>G4HEPEvtInterface</em></span> class object instantiated with a signal
|
|---|
| 229 | event file is accessed once. Thus, this event represents a typical
|
|---|
| 230 | signal event of LHC overlapping 20 minimum bias events. It should
|
|---|
| 231 | be noted that a simulation of event overlapping can be done by
|
|---|
| 232 | merging hits and/or digits associated with several events, and
|
|---|
| 233 | these events can be simulated independently. Digitization over
|
|---|
| 234 | multiple events will be mentioned in
|
|---|
| 235 | <a href="ch04s05.html" title="4.5.
|
|---|
| 236 | Digitization
|
|---|
| 237 | ">Section 4.5</a>.
|
|---|
| 238 | </p></div></div><div class="navfooter"><hr><table width="100%" summary="Navigation footer"><tr><td width="40%" align="left"><a accesskey="p" href="ch03s05.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><td width="20%" align="center"><a accesskey="u" href="ch03.html"><img src="AllResources/IconsGIF/up.gif" alt="Up"></a></td><td width="40%" align="right"> <a accesskey="n" href="ch03s07.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr><tr><td width="40%" align="left" valign="top">3.5.
|
|---|
| 239 | Event
|
|---|
| 240 | </td><td width="20%" align="center"><a accesskey="h" href="index.html"><img src="AllResources/IconsGIF/home.gif" alt="Home"></a></td><td width="40%" align="right" valign="top"> 3.7.
|
|---|
| 241 | Event Biasing Techniques
|
|---|
| 242 | </td></tr></table></div></body></html>
|
|---|