Geant4 11.2.2
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),
97 nSrcPairs(0),
101 deltasInside(false),
102 sigmasInside(false),
103 kaonsInside(false),
104 antikaonsInside(false),
105 lambdasInside(false),
106 forcedDeltasInside(false),
107 forcedDeltasOutside(false),
110 forcedSigmaOutside(false),
111 forcedStrangeInside(false),
112 emitLambda(0),
113 emitKaon(false),
114 clusterDecay(false),
122 nDecayAvatars(0),
125
126 {
127 std::fill_n(A, maxSizeParticles, 0);
128 std::fill_n(Z, maxSizeParticles, 0);
129 std::fill_n(S, maxSizeParticles, 0);
130 std::fill_n(J, maxSizeParticles, 0);
131 std::fill_n(PDGCode, maxSizeParticles, 0);
132 std::fill_n(ParticleBias, maxSizeParticles, (Float_t)0.0);
133 std::fill_n(EKin, maxSizeParticles, (Float_t)0.0);
134 std::fill_n(px, maxSizeParticles, (Float_t)0.0);
135 std::fill_n(py, maxSizeParticles, (Float_t)0.0);
136 std::fill_n(pz, maxSizeParticles, (Float_t)0.0);
137 std::fill_n(theta, maxSizeParticles, (Float_t)0.0);
138 std::fill_n(phi, maxSizeParticles, (Float_t)0.0);
139 std::fill_n(origin, maxSizeParticles, 0);
141 std::fill_n(parentResonanceID, maxSizeParticles, 0);
142 std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
143 std::fill_n(ARem, maxSizeRemnants, 0);
144 std::fill_n(ZRem, maxSizeRemnants, 0);
145 std::fill_n(SRem, maxSizeRemnants, 0);
146 std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
147 std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
148 std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
149 std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
150 std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
151 std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
152 std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
153 std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
154 std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
155 std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
156 std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
157 std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
158 std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
159 std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
160 }
161
162 /** \brief Number of the event */
164
165 /** \brief Maximum array size for remnants */
166 static const Short_t maxSizeRemnants = 10;
167
168 /** \brief Maximum array size for emitted particles */
169 static const Short_t maxSizeParticles = 1000;
170
171 /** \brief Number of particles in the final state */
173 /** \brief Sequential number of the event in the event loop */
175 /** \brief Particle mass number */
177 /** \brief Particle charge number */
179 /** \brief Particle strangeness number */
181 /** \brief Particle angular momemtum */
183 /** \brief PDG numbering of the particles */
185 /** \brief Event bias */
187 /** \brief Particle weight due to the bias */
189 /** \brief Particle kinetic energy [MeV] */
191 /** \brief Particle momentum, x component [MeV/c] */
193 /** \brief Particle momentum, y component [MeV/c] */
195 /** \brief Particle momentum, z component [MeV/c] */
197 /** \brief Particle momentum polar angle [radians] */
199 /** \brief Particle momentum azimuthal angle [radians] */
201 /** \brief Origin of the particle
202 *
203 * Should be -1 for cascade particles, or the number of the remnant for
204 * de-excitation particles. */
206 /** \brief Particle's parent resonance PDG code */
208 /** \brief Particle's parent resonance unique ID identifier */
210 /** \brief History of the particle
211 *
212 * Condensed information about the de-excitation chain of a particle. For
213 * cascade particles, it is just an empty string. For particles arising
214 * from the de-excitation of a cascade remnant, it is a string of
215 * characters. Each character represents one or more identical steps in
216 * the de-excitation process. The currently defined possible character
217 * values and their meanings are the following:
218 *
219 * e: evaporation product
220 * E: evaporation residue
221 * m: multifragmentation
222 * a: light partner in asymmetric fission or IMF emission
223 * A: heavy partner in asymmetric fission or IMF emission
224 * f: light partner in fission
225 * F: heavy partner in fission
226 * s: saddle-to-scission emission
227 * n: non-statistical emission (decay) */
228 std::vector<std::string> history;
229 /** \brief Number of remnants */
231 /** \brief Projectile particle type */
233 /** \brief Mass number of the target nucleus */
235 /** \brief Charge number of the target nucleus */
237 /** \brief Strangeness number of the target nucleus */
239 /** \brief Mass number of the projectile nucleus */
241 /** \brief Charge number of the projectile nucleus */
243 /** \brief Strangeness number of the projectile nucleus */
245 /** \brief Projectile kinetic energy given as input */
247 /** \brief Impact parameter [fm] */
249 /** \brief Number of accepted two-body collisions */
251 /** \brief Cascade stopping time [fm/c] */
253 /** \brief Energy-conservation balance [MeV] */
255 /** \brief First value for the energy-conservation balance [MeV] */
257 /** \brief Longitudinal momentum-conservation balance [MeV/c] */
259 /** \brief Transverse momentum-conservation balance [MeV/c] */
261 /** \brief Number of cascade particles */
263 /** \brief True if the event is transparent */
265 /** \brief True if annihilation at rest on a proton */
267 /** \brief True if annihilation at rest on a neutron */
269 /** \brief True if the event is a forced CN */
271 /** \brief True if the event is a nucleon absorption */
273 /** \brief True if the event is a pion absorption */
275 /** \brief Number of accepted Delta decays */
277 /** \brief Number of accepted SRC collisions */
279 /** \brief Number of src pairs */
281 /** \brief Number of two-body collisions blocked by Pauli or CDPP */
283 /** \brief Number of decays blocked by Pauli or CDPP */
285 /** \brief Effective (Coulomb-distorted) impact parameter [fm] */
287 /** \brief Event involved deltas in the nucleus at the end of the cascade */
289 /** \brief Event involved sigmas in the nucleus at the end of the cascade */
291 /** \brief Event involved kaons in the nucleus at the end of the cascade */
293 /** \brief Event involved antikaons in the nucleus at the end of the cascade */
295 /** \brief Event involved lambdas in the nucleus at the end of the cascade */
297 /** \brief Event involved forced delta decays inside the nucleus */
299 /** \brief Event involved forced delta decays outside the nucleus */
301 /** \brief Event involved forced eta/omega decays outside the nucleus */
303 /** \brief Event involved forced strange absorption inside the nucleus */
305 /** \brief Event involved forced Sigma Zero decays outside the nucleus */
307 /** \brief Event involved forced antiKaon/Sigma absorption inside the nucleus */
309 /** \brief Number of forced Lambda emit out of the nucleus */
311 /** \brief Event involved forced Kaon emission */
313 /** \brief Event involved cluster decay */
315 /** \brief Time of the first collision [fm/c] */
317 /** \brief Cross section of the first collision (mb) */
319 /** \brief Position of the spectator on the first collision (fm) */
321 /** \brief Momentum of the spectator on the first collision (fm) */
323 /** \brief True if the first collision was elastic */
325 /** \brief Number of reflection avatars */
327 /** \brief Number of collision avatars */
329 /** \brief Number of decay avatars */
331 /** \brief Number of dynamical spectators that were merged back into the projectile remnant */
333 /** \brief Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. */
335 /** \brief Emission time [fm/c] */
337 /** \brief Remnant mass number */
339 /** \brief Remnant charge number */
341 /** \brief Remnant strangeness number */
343 /** \brief Remnant excitation energy [MeV] */
345 /** \brief Remnant spin [\f$\hbar\f$] */
347 /** \brief Remnant kinetic energy [MeV] */
349 /** \brief Remnant momentum, x component [MeV/c] */
351 /** \brief Remnant momentum, y component [MeV/c] */
353 /** \brief Remnant momentum, z component [MeV/c] */
355 /** \brief Remnant momentum polar angle [radians] */
357 /** \brief Remnant momentum azimuthal angle [radians] */
359 /** \brief Remnant angular momentum, x component [\f$\hbar\f$] */
361 /** \brief Remnant angular momentum, y component [\f$\hbar\f$] */
363 /** \brief Remnant angular momentum, z component [\f$\hbar\f$] */
365 /** \brief Particle kinetic energy, in inverse kinematics [MeV] */
367 /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */
369 /** \brief Particle momentum polar angle, in inverse kinematics [radians] */
371
372 /** \brief Reset the EventInfo members */
373 void reset() {
374 nParticles = 0;
375 event = 0;
376 eventBias = (Float_t)0.0;
377 history.clear();
378 nRemnants = 0;
379 projectileType = 0;
380 At = 0;
381 Zt = 0;
382 St = 0;
383 Ap = 0;
384 Zp = 0;
385 Sp = 0;
386 Ep = (Float_t)0.0;
388 nCollisions = 0;
389 stoppingTime = (Float_t)0.0;
390 EBalance = (Float_t)0.0;
391 firstEBalance = (Float_t)0.0;
392 pLongBalance = (Float_t)0.0;
393 pTransBalance = (Float_t)0.0;
395 transparent = false;
396 annihilationP = false;
397 annihilationN = false;
398 forcedCompoundNucleus = false;
399 nucleonAbsorption = false;
400 pionAbsorption = false;
401 nDecays = 0;
402 nSrcCollisions = 0;
403 nSrcPairs = 0;
405 nBlockedDecays = 0;
407 deltasInside = false;
408 sigmasInside = false;
409 kaonsInside = false;
410 antikaonsInside = false;
411 lambdasInside = false;
412 forcedDeltasInside = false;
413 forcedDeltasOutside = false;
416 forcedSigmaOutside = false;
417 forcedStrangeInside = false;
418 emitLambda = 0;
419 emitKaon = false;
420 clusterDecay = false;
428 nDecayAvatars = 0;
431
432 }
433
434 /// \brief Move a remnant to the particle array
435 void remnantToParticle(const G4int remnantIndex);
436
437 /// \brief Fill the variables describing the reaction in inverse kinematics
438 void fillInverseKinematics(const Double_t gamma);
439 };
440}
441
442#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.
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.
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].
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].
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].
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 [ ].
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.
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].
Float_t firstCollisionTime
Time of the first collision [fm/c].
Short_t Zt
Charge number of the target nucleus.
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