source: trunk/documents/UserDoc/DocBookUsersGuides/ForApplicationDeveloper/xml/Fundamentals/eventGenerator.xml @ 1013

Last change on this file since 1013 was 904, checked in by garnier, 16 years ago

ajout de la doc

File size: 11.4 KB
Line 
1<!-- ******************************************************** -->
2<!--                                                          -->
3<!--  [History]                                               -->
4<!--    Converted to DocBook: Katsuya Amako, Aug-2006         -->
5<!--    Changed by: Katsuya Amako, 21-Sep-1998                -->
6<!--    Changed by: Dennis Wright, 25-Jun-2002                -->
7<!--    Proof read by: Joe Chuma,  28-Jun-1999                -->
8<!--                                                          -->
9<!-- ******************************************************** -->
10
11
12<!-- ******************* Section (Level#1) ****************** -->
13<sect1 id="sect.EventGen">
14<title>
15Event Generator Interface
16</title>
17
18<!-- ******************* Section (Level#2) ****************** -->
19<sect2 id="sect.EventGen.StructPriEvt">
20<title>
21Structure of a primary event
22</title>
23
24<!-- ******************* Section (Level#3) ****************** -->
25<sect3 id="sect.EventGen.StructPriEvt.VtxPart">
26<title>
27Primary vertex and primary particle
28</title>
29
30<para>
31The <emphasis>G4Event</emphasis> class object should have a set of primary
32particles when it is sent to <emphasis>G4EventManager</emphasis> via
33<literal>processOneEvent()</literal> method. It is the mandate of your
34<emphasis>G4VUserPrimaryGeneratorAction</emphasis> concrete class to send primary
35particles to the <emphasis>G4Event</emphasis> object.
36</para>
37
38<para>
39The <emphasis>G4PrimaryParticle</emphasis> class represents a primary particle
40with which Geant4 starts simulating an event. This class object has
41information on particle type and its three momenta. The positional
42and time information of primary particle(s) are stored in the
43<emphasis>G4PrimaryVertex</emphasis> class object and, thus, this class object
44can have one or more <emphasis>G4PrimaryParticle</emphasis> class objects which
45share the same vertex. As shown in Fig.?.?, primary vertexes and
46primary particles are associated with the <emphasis>G4Event</emphasis> object by
47a form of linked list.
48</para>
49
50<para>
51A concrete class of <emphasis>G4VPrimaryGenerator</emphasis>, the
52<emphasis>G4PrimaryParticle</emphasis> object is constructed with either a
53pointer to <emphasis>G4ParticleDefinition</emphasis> or an integer number which
54represents P.D.G. particle code. For the case of some artificial
55particles, e.g., geantino, optical photon, etc., or exotic nuclear
56fragments, which the P.D.G. particle code does not cover, the
57<emphasis>G4PrimaryParticle</emphasis> should be constructed by
58<emphasis>G4ParticleDefinition</emphasis> pointer. On the other hand, elementary
59particles with very short life time, e.g., weak bosons, or
60quarks/gluons, can be instantiated as <emphasis>G4PrimaryParticle</emphasis>
61objects using the P.D.G. particle code. It should be noted that,
62even though primary particles with such a very short life time are
63defined, Geant4 will simulate only the particles which are defined
64as <emphasis>G4ParticleDefinition</emphasis> class objects. Other primary
65particles will be simply ignored by <emphasis>G4EventManager</emphasis>. But it
66may still be useful to construct such "intermediate" particles for
67recording the origin of the primary event.
68</para>
69
70</sect3>
71
72<!-- ******************* Section (Level#3) ****************** -->
73<sect3 id="sect.EventGen.StructPriEvt.ForcedDecay">
74<title>
75Forced decay channel
76</title>
77
78<para>
79The <emphasis>G4PrimaryParticle</emphasis> class object can have a list of its
80daughter particles. If the parent particle is an "intermediate"
81particle, which Geant4 does not have a corresponding
82<emphasis>G4ParticleDefinition</emphasis>, this parent particle is ignored and
83daughters are assumed to start from the vertex with which their
84parent is associated. For example, a Z boson is associated with a
85vertex and it has positive and negative muons as its daughters,
86these muons will start from that vertex.
87</para>
88
89<para>
90There are some kinds of particles which should fly some
91reasonable distances and, thus, should be simulated by Geant4, but
92you still want to follow the decay channel generated by an event
93generator. A typical case of these particles is B meson. Even for
94the case of a primary particle which has a corresponding
95<emphasis>G4ParticleDefinition</emphasis>, it can have daughter primary
96particles. Geant4 will trace the parent particle until it comes to
97decay, obeying multiple scattering, ionization loss, rotation with
98the magnetic field, etc. according to its particle type. When the
99parent comes to decay, instead of randomly choosing its decay
100channel, it follows the "pre-assigned" decay channel. To conserve
101the energy and the momentum of the parent, daughters will be
102Lorentz transformed according to their parent's frame.
103</para>
104
105</sect3>
106</sect2>
107
108<!-- ******************* Section (Level#2) ****************** -->
109<sect2 id="sect.EventGen.InterPriGen">
110<title>
111Interface to a primary generator
112</title>
113
114<!-- ******************* Section (Level#3) ****************** -->
115<sect3 id="sect.EventGen.InterPriGen.HEPEvt">
116<title>
117G4HEPEvtInterface
118</title>
119
120<para>
121Unfortunately, almost all event generators presently in use,
122commonly are written in FORTRAN. For Geant4, it was decided to not
123link with any FORTRAN program or library, even though the C++
124language syntax itself allows such a link. Linking to a FORTRAN
125package might be convenient in some cases, but we will lose many
126advantages of object-oriented features of C++, such as robustness.
127Instead, Geant4 provides an ASCII file interface for such event
128generators.
129</para>
130
131<para>
132<emphasis>G4HEPEvtInterface</emphasis> is one of <emphasis>G4VPrimaryGenerator</emphasis>
133concrete class and thus it can be used in your
134<emphasis>G4VUserPrimaryGeneratorAction</emphasis> concrete class.
135<emphasis>G4HEPEvtInterface</emphasis> reads an ASCII file produced by an event
136generator and reproduces <emphasis>G4PrimaryParticle</emphasis> objects
137associated with a <emphasis>G4PrimaryVertex</emphasis> object. It reproduces a
138full production chain of the event generator, starting with primary
139quarks, etc. In other words, <emphasis>G4HEPEvtInterface</emphasis> converts
140information stored in the <literal>/HEPEVT/</literal> common block to an
141object-oriented data structure. Because the <literal>/HEPEVT/</literal>
142common block is commonly used by almost all event generators
143written in FORTRAN, <emphasis>G4HEPEvtInterface</emphasis> can interface to
144almost all event generators currently used in the HEP community.
145The constructor of <emphasis>G4HEPEvtInterface</emphasis> takes the file name.
146<xref linkend="programlist_EventGen_1" /> shows an example how to use
147<emphasis>G4HEPEvtInterface</emphasis>. Note that an event generator is not
148assumed to give a place of the primary particles, the interaction
149point must be set before invoking <emphasis>GeneratePrimaryVertex()</emphasis>
150method.
151</para>
152
153<example id="programlist_EventGen_1">
154<title>
155An example code for <emphasis>G4HEPEvtInterface</emphasis>
156</title>
157<programlisting>
158#ifndef ExN04PrimaryGeneratorAction_h
159#define ExN04PrimaryGeneratorAction_h 1
160
161#include "G4VUserPrimaryGeneratorAction.hh"
162#include "globals.hh"
163
164class G4VPrimaryGenerator;
165class G4Event;
166
167class ExN04PrimaryGeneratorAction : public G4VUserPrimaryGeneratorAction
168{
169  public:
170    ExN04PrimaryGeneratorAction();
171    ~ExN04PrimaryGeneratorAction();
172
173  public:
174    void GeneratePrimaries(G4Event* anEvent);
175
176  private:
177    G4VPrimaryGenerator* HEPEvt;
178};
179
180#endif
181
182
183#include "ExN04PrimaryGeneratorAction.hh"
184
185#include "G4Event.hh"
186#include "G4HEPEvtInterface.hh"
187
188ExN04PrimaryGeneratorAction::ExN04PrimaryGeneratorAction()
189{
190  HEPEvt = new G4HEPEvtInterface("pythia_event.data");
191}
192
193ExN04PrimaryGeneratorAction::~ExN04PrimaryGeneratorAction()
194{
195  delete HEPEvt;
196}
197
198void ExN04PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
199{
200  HEPEvt-&gt;SetParticlePosition(G4ThreeVector(0.*cm,0.*cm,0.*cm));
201  HEPEvt-&gt;GeneratePrimaryVertex(anEvent);
202}
203</programlisting>
204</example>
205
206</sect3>
207
208<!-- ******************* Section (Level#3) ****************** -->
209<sect3 id="sect.EventGen.InterPriGen.FormASCII">
210<title>
211Format of the ASCII file
212</title>
213
214<para>
215An ASCII file, which will be fed by <emphasis>G4HEPEvtInterface</emphasis>
216should have the following format.
217
218<itemizedlist spacing="compact">
219  <listitem><para>
220    The first line of each primary event should be an integer which
221    represents the number of the following lines of primary
222    particles.
223  </para></listitem>
224  <listitem><para>
225    Each line in an event corresponds to a particle in the
226    <literal>/HEPEVT/</literal> common. Each line has <literal>ISTHEP, IDHEP,
227    JDAHEP(1), JDAHEP(2), PHEP(1), PHEP(2), PHEP(3), PHEP(5).</literal>
228    Refer to the <literal>/HEPEVT/</literal> manual for the meanings of these
229    variables.
230  </para></listitem>
231</itemizedlist>
232</para>
233
234<para>
235<xref linkend="programlist_EventGen_2" /> shows an example FORTRAN code to
236generate an ASCII file.
237</para>
238
239<example id="programlist_EventGen_2">
240<title>
241A FORTRAN example using the <literal>/HEPEVT/</literal> common.
242</title>
243<programlisting>
244***********************************************************
245      SUBROUTINE HEP2G4
246*
247* Convert /HEPEVT/ event structure to an ASCII file
248* to be fed by G4HEPEvtInterface
249*
250***********************************************************
251      PARAMETER (NMXHEP=2000)
252      COMMON/HEPEVT/NEVHEP,NHEP,ISTHEP(NMXHEP),IDHEP(NMXHEP),
253     &gt;JMOHEP(2,NMXHEP),JDAHEP(2,NMXHEP),PHEP(5,NMXHEP),VHEP(4,NMXHEP)
254      DOUBLE PRECISION PHEP,VHEP
255*
256      WRITE(6,*) NHEP
257      DO IHEP=1,NHEP
258       WRITE(6,10)
259     &gt;  ISTHEP(IHEP),IDHEP(IHEP),JDAHEP(1,IHEP),JDAHEP(2,IHEP),
260     &gt;  PHEP(1,IHEP),PHEP(2,IHEP),PHEP(3,IHEP),PHEP(5,IHEP)
26110     FORMAT(4I10,4(1X,D15.8))
262      ENDDO
263*
264      RETURN
265      END
266</programlisting>
267</example>
268
269</sect3>
270
271<!-- ******************* Section (Level#3) ****************** -->
272<sect3 id="sect.EventGen.InterPriGen.FutureInt">
273<title>
274Future interface to the new generation generators
275</title>
276
277<para>
278Several activities have already been started for developing
279object-oriented event generators. Such new generators can be easily
280linked and used with a Geant4 based simulation. Furthermore, we
281need not distinguish a primary generator from the physics processes
282used in Geant4. Future generators can be a kind of physics process
283plugged-in by inheriting <emphasis>G4VProcess</emphasis>.
284</para>
285
286</sect3>
287</sect2>
288
289
290<!-- ******************* Section (Level#2) ****************** -->
291<sect2 id="sect.EventGen.EvtOverlap">
292<title>
293Event overlap using multiple generators
294</title>
295
296<para>
297Your <emphasis>G4VUserPrimaryGeneratorAction</emphasis> concrete class can have
298more than one <emphasis>G4VPrimaryGenerator</emphasis> concrete class. Each
299<emphasis>G4VPrimaryGenerator</emphasis> concrete class can be accessed more than
300once per event. Using these class objects, one event can have more
301than one primary event.
302</para>
303
304<para>
305One possible use is the following. Within an event, a
306<emphasis>G4HEPEvtInterface</emphasis> class object instantiated with a minimum
307bias event file is accessed 20 times and another
308<emphasis>G4HEPEvtInterface</emphasis> class object instantiated with a signal
309event file is accessed once. Thus, this event represents a typical
310signal event of LHC overlapping 20 minimum bias events. It should
311be noted that a simulation of event overlapping can be done by
312merging hits and/or digits associated with several events, and
313these events can be simulated independently. Digitization over
314multiple events will be mentioned in
315<xref linkend="sect.Digi" />.
316</para>
317
318
319</sect2>
320</sect1>
Note: See TracBrowser for help on using the repository browser.