Geant4 9.6.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// Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28// Davide Mancusi, CEA
29// Alain Boudard, CEA
30// Sylvie Leray, CEA
31// Joseph Cugnon, University of Liege
32//
33// INCL++ revision: v5.1.8
34//
35#define INCLXX_IN_GEANT4_MODE 1
36
37#include "globals.hh"
38
39/** \file G4INCLEventInfo.hh
40 * \brief Simple container for output of event results.
41 *
42 * Contains the results of an INCL cascade.
43 *
44 * \date 21 January 2011
45 * \author Davide Mancusi
46 */
47
48#ifndef G4INCLEVENTINFO_HH
49#define G4INCLEVENTINFO_HH 1
50
51#include "G4INCLParticleType.hh"
52#ifdef INCL_ROOT_USE
53#include <Rtypes.h>
54#endif
55#include <string>
56#include <vector>
57#include <algorithm>
58
59namespace G4INCL {
60#ifndef INCL_ROOT_USE
61 typedef G4int Int_t;
62 typedef short Short_t;
65 typedef G4bool Bool_t;
66#endif
67
68 struct EventInfo {
71 At(0), Zt(0), Ap(0), Zp(0),
72 Ep(0.),
74 EBalance(0.0), pLongBalance(0.0), pTransBalance(0.0),
76 transparent(true),
78 nucleonAbsorption(false), pionAbsorption(false), nDecays(0),
81 deltasInside(false),
82 forcedDeltasInside(false),
84 clusterDecay(false),
91 {
92 std::fill_n(ARem, maxSizeRemnants, 0);
93 std::fill_n(ZRem, maxSizeRemnants, 0);
94 std::fill_n(EStarRem, maxSizeRemnants, ((Float_t)0.));
95 std::fill_n(JRem, maxSizeRemnants, ((Float_t)0.));
96 std::fill_n(EKinRem, maxSizeRemnants, ((Float_t)0.));
97 std::fill_n(pxRem, maxSizeRemnants, ((Float_t)0.));
98 std::fill_n(pyRem, maxSizeRemnants, ((Float_t)0.));
99 std::fill_n(pzRem, maxSizeRemnants, ((Float_t)0.));
100 std::fill_n(thetaRem, maxSizeRemnants, ((Float_t)0.));
101 std::fill_n(phiRem, maxSizeRemnants, ((Float_t)0.));
102 std::fill_n(jxRem, maxSizeRemnants, ((Float_t)0.));
103 std::fill_n(jyRem, maxSizeRemnants, ((Float_t)0.));
104 std::fill_n(jzRem, maxSizeRemnants, ((Float_t)0.));
105
106 std::fill_n(A, maxSizeParticles, 0);
107 std::fill_n(Z, maxSizeParticles, 0);
108 std::fill_n(emissionTime, maxSizeParticles, ((Float_t)0.));
109 std::fill_n(EKin, maxSizeParticles, ((Float_t)0.));
110 std::fill_n(px, maxSizeParticles, ((Float_t)0.));
111 std::fill_n(py, maxSizeParticles, ((Float_t)0.));
112 std::fill_n(pz, maxSizeParticles, ((Float_t)0.));
113 std::fill_n(theta, maxSizeParticles, ((Float_t)0.));
114 std::fill_n(phi, maxSizeParticles, ((Float_t)0.));
115 std::fill_n(origin, maxSizeParticles, 0);
116 };
117
118 /** \brief Number of the event */
120 /** \brief Protjectile particle type */
122
123 /** \brief Mass number of the target nucleus */
125 /** \brief Charge number of the target nucleus */
127
128 /** \brief Mass number of the projectile nucleus */
130 /** \brief Charge number of the projectile nucleus */
132 /** \brief Projectile kinetic energy given as input */
134
135 /** \brief Impact parameter [fm] */
137 /** \brief Number of accepted two-body collisions */
139 /** \brief Cascade stopping time [fm/c] */
141
142 /** \brief Energy-conservation balance [MeV] */
144 /** \brief Longitudinal momentum-conservation balance [MeV/c] */
146 /** \brief Transverse momentum-conservation balance [MeV/c] */
148
149 /** \brief Number of cascade particles */
151 /** \brief Number of remnants */
153 /** \brief Total number of emitted particles */
155
156 /** \brief True if the event is transparent */
158 /** \brief True if the event is a forced CN */
160 /** \brief True if the event is absorption */
162 /** \brief True if the event is absorption */
164 /** \brief Number of accepted Delta decays */
166 /** \brief Number of two-body collisions blocked by Pauli or CDPP */
168 /** \brief Number of decays blocked by Pauli or CDPP */
170 /** \brief Number of reflection avatars */
171 /** \brief Effective (Coulomb-distorted) impact parameter [fm] */
173
174 /// \brief Event involved deltas in the nucleus at the end of the cascade
176 /// \brief Event involved forced delta decays inside the nucleus
178 /// \brief Event involved forced delta decays outside the nucleus
180 /// \brief Event involved cluster decay
182
183 /** \brief Time of the first collision [fm/c] */
185 /** \brief Cross section of the first collision (mb) */
187
189 /** \brief Number of collision avatars */
191 /** \brief Number of decay avatars */
193
194 /// \brief Number of dynamical spectators that were merged back into the projectile remnant
196
197 /** \brief Maximum array size for remnants */
198 static const Short_t maxSizeRemnants = 10;
199 /** \brief Remnant mass number */
201 /** \brief Remnant charge number */
203 /** \brief Remnant excitation energy [MeV] */
205 /** \brief Remnant spin [\f$\hbar\f$] */
207 /** \brief Remnant kinetic energy [MeV] */
209 /** \brief Remnant momentum, x component [MeV/c] */
211 /** \brief Remnant momentum, y component [MeV/c] */
213 /** \brief Remnant momentum, z component [MeV/c] */
215 /** \brief Remnant momentum polar angle [radians] */
217 /** \brief Remnant momentum azimuthal angle [radians] */
219 /** \brief Remnant angular momentum, x component [hbar] */
221 /** \brief Remnant angular momentum, y component [hbar] */
223 /** \brief Remnant angular momentum, z component [hbar] */
225
226 /** \brief Maximum array size for emitted particles */
227 static const Short_t maxSizeParticles = 1000;
228 /** \brief Particle mass number */
230 /** \brief Particle charge number */
232 /** \brief Emission time [fm/c] */
234 /** \brief Particle kinetic energy [MeV] */
236 /** \brief Particle momentum, x component [MeV/c] */
238 /** \brief Particle momentum, y component [MeV/c] */
240 /** \brief Particle momentum, z component [MeV/c] */
242 /** \brief Particle momentum polar angle [radians] */
244 /** \brief Particle momentum azimuthal angle [radians] */
246 /** \brief Origin of the particle
247 *
248 * Should be -1 for cascade particles, or the number of the remnant for
249 * de-excitation particles.
250 *
251 */
253 /** \brief History of the particle
254 *
255 * Condensed information about the de-excitation chain of a particle. For
256 * cascade particles, it is just an empty string. For particles arising
257 * from the de-excitation of a cascade remnant, it is a string of
258 * characters. Each character represents one or more identical steps in
259 * the de-excitation process. The currently defined possible character
260 * values and their meanings are the following:
261 *
262 * e: evaporation product
263 * E: evaporation residue
264 * m: multifragmentation
265 * a: light partner in asymmetric fission or IMF emission
266 * A: heavy partner in asymmetric fission or IMF emission
267 * f: light partner in fission
268 * F: heavy partner in fission
269 * s: saddle-to-scission emission
270 * n: non-statistical emission (decay)
271 */
272 std::vector<std::string> history;
273
274#ifdef INCL_INVERSE_KINEMATICS
275 /** \brief Particle kinetic energy, in inverse kinematics [MeV] */
276 Float_t EKinPrime[maxSizeParticles];
277 /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */
278 Float_t pzPrime[maxSizeParticles];
279 /** \brief Particle momentum polar angle, in inverse kinematics [radians] */
280 Float_t thetaPrime[maxSizeParticles];
281#endif // INCL_INVERSE_KINEMATICS
282
283 /** \brief Reset the EventInfo members */
284 void reset() {
285 Ap = 0;
286 Zp = 0;
287 At = 0;
288 Zt = 0;
289 impactParameter = 0.0;
291 stoppingTime = 0.0;
292 EBalance = 0.0;
293 pLongBalance = 0.0;
294 pTransBalance = 0.0;
295 nCollisions = 0;
297 nDecays = 0;
299 nDecays = 0;
301 nRemnants = 0;
302 nParticles = 0;
303 transparent = true;
304 forcedCompoundNucleus = false;
305 nucleonAbsorption = false;
306 pionAbsorption = false;
307 forcedDeltasInside = false;
308 forcedDeltasOutside = false;
309 deltasInside = false;
310 clusterDecay = false;
312 }
313
314#ifdef INCL_INVERSE_KINEMATICS
315 void fillInverseKinematics(const Double_t gamma);
316#endif // INCL_INVERSE_KINEMATICS
317 };
318}
319
320#endif /* G4INCLEVENTINFO_HH */
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
short Short_t
G4bool Bool_t
G4float Float_t
G4int Int_t
G4double Double_t
Bool_t pionAbsorption
True if the event is absorption.
void reset()
Reset the EventInfo members.
Int_t nCollisionAvatars
Number of collision avatars.
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 nDecays
Number of accepted Delta decays.
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
static Int_t eventNumber
Number of the event.
Float_t Ep
Projectile kinetic energy given as input.
Short_t nCascadeParticles
Number of cascade particles.
Bool_t transparent
True if the event is transparent.
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Short_t A[maxSizeParticles]
Particle mass number.
Bool_t deltasInside
Event involved deltas in the nucleus at the end of the cascade.
Float_t EBalance
Energy-conservation balance [MeV].
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [hbar].
Bool_t nucleonAbsorption
True if the event is absorption.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
static const Short_t maxSizeRemnants
Maximum array size for remnants.
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
Float_t effectiveImpactParameter
Number of reflection avatars.
Float_t stoppingTime
Cascade stopping time [fm/c].
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
Short_t Ap
Mass number of the projectile nucleus.
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.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [hbar].
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
Bool_t clusterDecay
Event involved cluster decay.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t firstCollisionTime
Time of the first collision [fm/c].
Short_t Zt
Charge number of the target nucleus.
Int_t nParticles
Total number of emitted particles.
Float_t phiRem[maxSizeRemnants]
Remnant momentum azimuthal angle [radians].
ParticleType projectileType
Protjectile particle type.
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
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.
Int_t nRemnants
Number of remnants.
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant.
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 [hbar].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
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.