Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLEventInfo.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// INCL++ intra-nuclear cascade model
27// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
38/** \file G4INCLEventInfo.hh
39 * \brief Simple container for output of event results.
40 *
41 * Contains the results of an INCL cascade.
42 *
43 * \date 21 January 2011
44 * \author Davide Mancusi
45 */
46
47#ifndef G4INCLEVENTINFO_HH_HH
48#define G4INCLEVENTINFO_HH_HH 1
49
50#include "G4INCLParticleType.hh"
51#ifdef INCL_ROOT_USE
52#include <Rtypes.h>
53#endif
54#include <string>
55#include <vector>
56#include <algorithm>
57
58namespace G4INCL {
59#ifndef INCL_ROOT_USE
60 typedef G4int Int_t;
61 typedef short Short_t;
64 typedef G4bool Bool_t;
65#endif
66
67 struct EventInfo {
69 nParticles(0),
70 event(0),
71 eventBias((Float_t)0.0),
72 nRemnants(0),
74 At(0),
75 Zt(0),
76 St(0),
77 Ap(0),
78 Zp(0),
79 Sp(0),
80 Ep((Float_t)0.0),
82 nCollisions(0),
84 EBalance((Float_t)0.0),
89 transparent(false),
90 annihilationP(false),
91 annihilationN(false),
93 nucleonAbsorption(false),
94 pionAbsorption(false),
95 nDecays(0),
96 fission(false),
97 fissmode(0),
98 EStarFis((Float_t)0.0),
99 ASad(0),
100 ZSad(0),
102 nSrcPairs(0),
106 deltasInside(false),
107 sigmasInside(false),
108 kaonsInside(false),
109 antinucleonsInside(false),
110 antikaonsInside(false),
111 lambdasInside(false),
112 forcedDeltasInside(false),
113 forcedDeltasOutside(false),
116 forcedSigmaOutside(false),
117 forcedStrangeInside(false),
118 emitLambda(0),
120 emitKaon(false),
121 emitAntinucleon(false),
122 clusterDecay(false),
130 nDecayAvatars(0),
133
134 {
135 std::fill_n(A, maxSizeParticles, 0);
136 std::fill_n(Z, maxSizeParticles, 0);
137 std::fill_n(S, maxSizeParticles, 0);
138 std::fill_n(J, maxSizeParticles, 0);
139 std::fill_n(PDGCode, maxSizeParticles, 0);
140 std::fill_n(ParticleBias, maxSizeParticles, (Float_t)0.0);
141 std::fill_n(EKin, maxSizeParticles, (Float_t)0.0);
142 std::fill_n(px, maxSizeParticles, (Float_t)0.0);
143 std::fill_n(py, maxSizeParticles, (Float_t)0.0);
144 std::fill_n(pz, maxSizeParticles, (Float_t)0.0);
145 std::fill_n(theta, maxSizeParticles, (Float_t)0.0);
146 std::fill_n(phi, maxSizeParticles, (Float_t)0.0);
147 std::fill_n(origin, maxSizeParticles, 0);
149 std::fill_n(parentResonanceID, maxSizeParticles, 0);
150 std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
151 std::fill_n(ARem, maxSizeRemnants, 0);
152 std::fill_n(ZRem, maxSizeRemnants, 0);
153 std::fill_n(SRem, maxSizeRemnants, 0);
154 std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
155 std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
156 std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
157 std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
158 std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
159 std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
160 std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
161 std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
162 std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
163 std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
164 std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
165 std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
166 std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
167 std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
168 }
169
170 /** \brief Number of the event */
172
173 /** \brief Maximum array size for remnants */
174 static const Short_t maxSizeRemnants = 10;
175
176 /** \brief Maximum array size for emitted particles */
177 static const Short_t maxSizeParticles = 1000;
178
179 /** \brief Number of particles in the final state */
181 /** \brief Sequential number of the event in the event loop */
183 /** \brief Particle mass number */
185 /** \brief Particle charge number */
187 /** \brief Particle strangeness number */
189 /** \brief Particle angular momemtum */
191 /** \brief PDG numbering of the particles */
193 /** \brief Event bias */
195 /** \brief Particle weight due to the bias */
197 /** \brief Particle kinetic energy [MeV] */
199 /** \brief Particle momentum, x component [MeV/c] */
201 /** \brief Particle momentum, y component [MeV/c] */
203 /** \brief Particle momentum, z component [MeV/c] */
205 /** \brief Particle momentum polar angle [radians] */
207 /** \brief Particle momentum azimuthal angle [radians] */
209 /** \brief Origin of the particle
210 *
211 * Should be -1 for cascade particles, or the number of the remnant for
212 * de-excitation particles. */
214 /** \brief Particle's parent resonance PDG code */
216 /** \brief Particle's parent resonance unique ID identifier */
218 /** \brief History of the particle
219 *
220 * Condensed information about the de-excitation chain of a particle. For
221 * cascade particles, it is just an empty string. For particles arising
222 * from the de-excitation of a cascade remnant, it is a string of
223 * characters. Each character represents one or more identical steps in
224 * the de-excitation process. The currently defined possible character
225 * values and their meanings are the following:
226 *
227 * e: evaporation product
228 * E: evaporation residue
229 * m: multifragmentation
230 * a: light partner in asymmetric fission or IMF emission
231 * A: heavy partner in asymmetric fission or IMF emission
232 * f: light partner in fission
233 * F: heavy partner in fission
234 * s: saddle-to-scission emission
235 * n: non-statistical emission (decay) */
236 std::vector<std::string> history;
237 /** \brief Number of remnants */
239 /** \brief Projectile particle type */
241 /** \brief Mass number of the target nucleus */
243 /** \brief Charge number of the target nucleus */
245 /** \brief Strangeness number of the target nucleus */
247 /** \brief Mass number of the projectile nucleus */
249 /** \brief Charge number of the projectile nucleus */
251 /** \brief Strangeness number of the projectile nucleus */
253 /** \brief Projectile kinetic energy given as input */
255 /** \brief Impact parameter [fm] */
257 /** \brief Number of accepted two-body collisions */
259 /** \brief Cascade stopping time [fm/c] */
261 /** \brief Energy-conservation balance [MeV] */
263 /** \brief First value for the energy-conservation balance [MeV] */
265 /** \brief Longitudinal momentum-conservation balance [MeV/c] */
267 /** \brief Transverse momentum-conservation balance [MeV/c] */
269 /** \brief Number of cascade particles */
271 /** \brief True if the event is transparent */
273 /** \brief True if annihilation at rest on a proton */
275 /** \brief True if annihilation at rest on a neutron */
277 /** \brief True if the event is a forced CN */
279 /** \brief True if the event is a nucleon absorption */
281 /** \brief True if the event is a pion absorption */
283 /** \brief Number of accepted Delta decays */
285 /** \brief True if the event is fission */
287 /** \brief Fission mode */
289 /** \brief Excitation energy above fission barrier [MeV] */
291 /** \brief Mass number at saddle */
293 /** \brief Charge number at saddle */
295 /** \brief Mass number at scission */
296 std::vector<Int_t> ASci;
297 /** \brief Charge number at scission */
298 std::vector<Int_t> ZSci;
299 /** \brief Number of accepted SRC collisions */
301 /** \brief Number of src pairs */
303 /** \brief Number of two-body collisions blocked by Pauli or CDPP */
305 /** \brief Number of decays blocked by Pauli or CDPP */
307 /** \brief Effective (Coulomb-distorted) impact parameter [fm] */
309 /** \brief Event involved deltas in the nucleus at the end of the cascade */
311 /** \brief Event involved sigmas in the nucleus at the end of the cascade */
313 /** \brief Event involved kaons in the nucleus at the end of the cascade */
315 /** \brief Event involved antinucleons in the nucleus at the end of the cascade */
317 /** \brief Event involved antikaons in the nucleus at the end of the cascade */
319 /** \brief Event involved lambdas in the nucleus at the end of the cascade */
321 /** \brief Event involved forced delta decays inside the nucleus */
323 /** \brief Event involved forced delta decays outside the nucleus */
325 /** \brief Event involved forced eta/omega decays outside the nucleus */
327 /** \brief Event involved forced strange absorption inside the nucleus */
329 /** \brief Event involved forced Sigma Zero decays outside the nucleus */
331 /** \brief Event involved forced antiKaon/Sigma absorption inside the nucleus */
333 /** \brief Number of forced Lambda emit out of the nucleus */
335 /** \brief Number of forced Antilambda emit out of the nucleus */
337 /** \brief Event involved forced Kaon emission */
339 /** \brief Event involved forced Antinucleon emission */
341 /** \brief Event involved cluster decay */
343 /** \brief Time of the first collision [fm/c] */
345 /** \brief Cross section of the first collision (mb) */
347 /** \brief Position of the spectator on the first collision (fm) */
349 /** \brief Momentum of the spectator on the first collision (fm) */
351 /** \brief True if the first collision was elastic */
353 /** \brief Number of reflection avatars */
355 /** \brief Number of collision avatars */
357 /** \brief Number of decay avatars */
359 /** \brief Number of dynamical spectators that were merged back into the projectile remnant */
361 /** \brief Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. */
363 /** \brief Emission time [fm/c] */
365 /** \brief Remnant mass number */
367 /** \brief Remnant charge number */
369 /** \brief Remnant strangeness number */
371 /** \brief Remnant excitation energy [MeV] */
373 /** \brief Remnant spin [\f$\hbar\f$] */
375 /** \brief Remnant kinetic energy [MeV] */
377 /** \brief Remnant momentum, x component [MeV/c] */
379 /** \brief Remnant momentum, y component [MeV/c] */
381 /** \brief Remnant momentum, z component [MeV/c] */
383 /** \brief Remnant momentum polar angle [radians] */
385 /** \brief Remnant momentum azimuthal angle [radians] */
387 /** \brief Remnant angular momentum, x component [\f$\hbar\f$] */
389 /** \brief Remnant angular momentum, y component [\f$\hbar\f$] */
391 /** \brief Remnant angular momentum, z component [\f$\hbar\f$] */
393 /** \brief Particle kinetic energy, in inverse kinematics [MeV] */
395 /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */
397 /** \brief Particle momentum polar angle, in inverse kinematics [radians] */
399
400 /** \brief Reset the EventInfo members */
401 void reset() {
402 nParticles = 0;
403 event = 0;
404 eventBias = (Float_t)0.0;
405 history.clear();
406 nRemnants = 0;
407 projectileType = 0;
408 At = 0;
409 Zt = 0;
410 St = 0;
411 Ap = 0;
412 Zp = 0;
413 Sp = 0;
414 Ep = (Float_t)0.0;
416 nCollisions = 0;
417 stoppingTime = (Float_t)0.0;
418 EBalance = (Float_t)0.0;
419 firstEBalance = (Float_t)0.0;
420 pLongBalance = (Float_t)0.0;
421 pTransBalance = (Float_t)0.0;
423 transparent = false;
424 annihilationP = false;
425 annihilationN = false;
426 forcedCompoundNucleus = false;
427 nucleonAbsorption = false;
428 pionAbsorption = false;
429 nDecays = 0;
430 fission = false;
431 fissmode = 0;
432 EStarFis = (Float_t)0.0;
433 ASad = 0;
434 ZSad = 0;
435 ASci.clear();
436 ZSci.clear();
437 nSrcCollisions = 0;
438 nSrcPairs = 0;
440 nBlockedDecays = 0;
442 deltasInside = false;
443 sigmasInside = false;
444 kaonsInside = false;
445 antinucleonsInside = false;
446 antikaonsInside = false;
447 lambdasInside = false;
448 forcedDeltasInside = false;
449 forcedDeltasOutside = false;
452 forcedSigmaOutside = false;
453 forcedStrangeInside = false;
454 emitLambda = 0;
455 emitAntilambda = 0;
456 emitKaon = false;
457 emitAntinucleon = false;
458 clusterDecay = false;
466 nDecayAvatars = 0;
469
470 }
471
472 /// \brief Move a remnant to the particle array
473 void remnantToParticle(const G4int remnantIndex);
474
475 /// \brief Fill the variables describing the reaction in inverse kinematics
476 void fillInverseKinematics(const Double_t gamma);
477 };
478}
479
480#endif /* G4INCLEVENTINFO_HH_HH */
float G4float
Definition G4Types.hh:84
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
short Short_t
G4bool Bool_t
G4float Float_t
G4double Double_t
Bool_t pionAbsorption
True if the event is a pion absorption.
void reset()
Reset the EventInfo members.
Short_t S[maxSizeParticles]
Particle strangeness number.
Int_t nCollisionAvatars
Number of collision avatars.
Bool_t lambdasInside
Event involved lambdas in the nucleus at the end of the cascade.
Short_t origin[maxSizeParticles]
Origin of the particle.
Bool_t forcedDeltasOutside
Event involved forced delta decays outside the nucleus.
Short_t At
Mass number of the target nucleus.
Short_t Zp
Charge number of the projectile nucleus.
Int_t parentResonancePDGCode[maxSizeParticles]
Particle's parent resonance PDG code.
Int_t nDecays
Number of accepted Delta decays.
Bool_t annihilationP
True if annihilation at rest on a proton.
Float_t EStarFis
Excitation energy above fission barrier [MeV].
Int_t projectileType
Projectile particle type.
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
Float_t Ep
Projectile kinetic energy given as input.
Short_t nCascadeParticles
Number of cascade particles.
Short_t nRemnants
Number of remnants.
Float_t eventBias
Event bias.
Bool_t transparent
True if the event is transparent.
std::vector< Int_t > ZSci
Charge number at scission.
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Short_t A[maxSizeParticles]
Particle mass number.
Bool_t sigmasInside
Event involved sigmas in the nucleus at the end of the cascade.
Bool_t forcedSigmaOutside
Event involved forced Sigma Zero decays outside the nucleus.
Float_t firstEBalance
First value for the energy-conservation balance [MeV].
Bool_t deltasInside
Event involved deltas in the nucleus at the end of the cascade.
Float_t EBalance
Energy-conservation balance [MeV].
Short_t ASad
Mass number at saddle.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm).
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
Int_t parentResonanceID[maxSizeParticles]
Particle's parent resonance unique ID identifier.
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
Bool_t forcedPionResonancesOutside
Event involved forced eta/omega decays outside the nucleus.
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Bool_t antinucleonsInside
Event involved antinucleons in the nucleus at the end of the cascade.
Int_t emitLambda
Number of forced Lambda emit out of the nucleus.
static const Short_t maxSizeRemnants
Maximum array size for remnants.
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
std::vector< Int_t > ASci
Mass number at scission.
Short_t St
Strangeness number of the target nucleus.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
Bool_t kaonsInside
Event involved kaons in the nucleus at the end of the cascade.
Float_t stoppingTime
Cascade stopping time [fm/c].
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
Int_t nSrcPairs
Number of src pairs.
Short_t Ap
Mass number of the projectile nucleus.
Int_t nReflectionAvatars
Number of reflection avatars.
Int_t nSrcCollisions
Number of accepted SRC collisions.
Float_t firstCollisionSpectatorPosition
Position of the spectator on the first collision (fm).
Short_t Z[maxSizeParticles]
Particle charge number.
Bool_t forcedCompoundNucleus
True if the event is a forced CN.
Float_t JRem[maxSizeRemnants]
Remnant spin [ ].
Float_t thetaRem[maxSizeRemnants]
Remnant momentum polar angle [radians].
std::vector< std::string > history
History of the particle.
Short_t Sp
Strangeness number of the projectile nucleus.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Bool_t emitKaon
Event involved forced Kaon emission.
Short_t nParticles
Number of particles in the final state.
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
Bool_t fission
True if the event is fission.
Int_t event
Sequential number of the event in the event loop.
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
Short_t J[maxSizeParticles]
Particle angular momemtum.
Bool_t antikaonsInside
Event involved antikaons in the nucleus at the end of the cascade.
Bool_t clusterDecay
Event involved cluster decay.
Short_t ZSad
Charge number at saddle.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Float_t ParticleBias[maxSizeParticles]
Particle weight due to the bias.
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Short_t SRem[maxSizeRemnants]
Remnant strangeness number.
Bool_t annihilationN
True if annihilation at rest on a neutron.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
Short_t fissmode
Fission mode.
Float_t firstCollisionTime
Time of the first collision [fm/c].
Short_t Zt
Charge number of the target nucleus.
Int_t emitAntilambda
Number of forced Antilambda emit out of the nucleus.
Bool_t emitAntinucleon
Event involved forced Antinucleon emission.
Float_t phiRem[maxSizeRemnants]
Remnant momentum azimuthal angle [radians].
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
static G4ThreadLocal Int_t eventNumber
Number of the event.
Float_t impactParameter
Impact parameter [fm].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Int_t nCollisions
Number of accepted two-body collisions.
Bool_t absorbedStrangeParticle
Event involved forced strange absorption inside the nucleus.
Int_t PDGCode[maxSizeParticles]
PDG numbering of the particles.
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant.
void fillInverseKinematics(const Double_t gamma)
Fill the variables describing the reaction in inverse kinematics.
Int_t nBlockedCollisions
Number of two-body collisions blocked by Pauli or CDPP.
static const Short_t maxSizeParticles
Maximum array size for emitted particles.
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
void remnantToParticle(const G4int remnantIndex)
Move a remnant to the particle array.
Bool_t forcedStrangeInside
Event involved forced antiKaon/Sigma absorption inside the nucleus.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
Float_t firstCollisionXSec
Cross section of the first collision (mb).
Float_t pLongBalance
Longitudinal momentum-conservation balance [MeV/c].
Int_t nDecayAvatars
Number of decay avatars.
#define G4ThreadLocal
Definition tls.hh:77