Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 eventBias((Float_t)0.0),
71 nRemnants(0),
73 At(0),
74 Zt(0),
75 St(0),
76 Ap(0),
77 Zp(0),
78 Sp(0),
79 Ep((Float_t)0.0),
81 nCollisions(0),
83 EBalance((Float_t)0.0),
87 transparent(false),
89 nucleonAbsorption(false),
90 pionAbsorption(false),
91 nDecays(0),
95 deltasInside(false),
96 sigmasInside(false),
97 kaonsInside(false),
98 antikaonsInside(false),
99 lambdasInside(false),
100 forcedDeltasInside(false),
101 forcedDeltasOutside(false),
104 forcedSigmaOutside(false),
105 forcedStrangeInside(false),
106 emitLambda(0),
107 emitKaon(false),
108 clusterDecay(false),
116 nDecayAvatars(0),
119 event(0)
120
121 {
122 std::fill_n(A, maxSizeParticles, 0);
123 std::fill_n(Z, maxSizeParticles, 0);
124 std::fill_n(S, maxSizeParticles, 0);
125 std::fill_n(PDGCode, maxSizeParticles, 0);
126 std::fill_n(ParticleBias, maxSizeParticles, (Float_t)0.0);
127 std::fill_n(EKin, maxSizeParticles, (Float_t)0.0);
128 std::fill_n(px, maxSizeParticles, (Float_t)0.0);
129 std::fill_n(py, maxSizeParticles, (Float_t)0.0);
130 std::fill_n(pz, maxSizeParticles, (Float_t)0.0);
131 std::fill_n(theta, maxSizeParticles, (Float_t)0.0);
132 std::fill_n(phi, maxSizeParticles, (Float_t)0.0);
133 std::fill_n(origin, maxSizeParticles, 0);
135 std::fill_n(parentResonanceID, maxSizeParticles, 0);
136 std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
137 std::fill_n(ARem, maxSizeRemnants, 0);
138 std::fill_n(ZRem, maxSizeRemnants, 0);
139 std::fill_n(SRem, maxSizeRemnants, 0);
140 std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
141 std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
142 std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
143 std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
144 std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
145 std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
146 std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
147 std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
148 std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
149 std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
150 std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
151 std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
152 std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
153 std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
154 }
155
156 /** \brief Number of the event */
158
159 /** \brief Maximum array size for remnants */
160 static const Short_t maxSizeRemnants = 10;
161
162 /** \brief Maximum array size for emitted particles */
163 static const Short_t maxSizeParticles = 1000;
164
165 /** \brief Number of particles in the final state */
167 /** \brief Particle mass number */
169 /** \brief Particle charge number */
171 /** \brief Particle strangeness number */
173 /** \brief PDG numbering of the particles */
175 /** \brief Particle weight due to the bias */
177 /** \brief Event bias */
179 /** \brief Particle kinetic energy [MeV] */
181 /** \brief Particle momentum, x component [MeV/c] */
183 /** \brief Particle momentum, y component [MeV/c] */
185 /** \brief Particle momentum, z component [MeV/c] */
187 /** \brief Particle momentum polar angle [radians] */
189 /** \brief Particle momentum azimuthal angle [radians] */
191 /** \brief Origin of the particle
192 *
193 * Should be -1 for cascade particles, or the number of the remnant for
194 * de-excitation particles. */
196 /** \brief Particle's parent resonance PDG code */
198 /** \brief Particle's parent resonance unique ID identifier */
200 /** \brief Emission time [fm/c] */
202 /** \brief History of the particle
203 *
204 * Condensed information about the de-excitation chain of a particle. For
205 * cascade particles, it is just an empty string. For particles arising
206 * from the de-excitation of a cascade remnant, it is a string of
207 * characters. Each character represents one or more identical steps in
208 * the de-excitation process. The currently defined possible character
209 * values and their meanings are the following:
210 *
211 * e: evaporation product
212 * E: evaporation residue
213 * m: multifragmentation
214 * a: light partner in asymmetric fission or IMF emission
215 * A: heavy partner in asymmetric fission or IMF emission
216 * f: light partner in fission
217 * F: heavy partner in fission
218 * s: saddle-to-scission emission
219 * n: non-statistical emission (decay) */
220 std::vector<std::string> history;
221 /** \brief Number of remnants */
223 /** \brief Remnant mass number */
225 /** \brief Remnant charge number */
227 /** \brief Remnant strangeness number */
229 /** \brief Remnant excitation energy [MeV] */
231 /** \brief Remnant spin [\f$\hbar\f$] */
233 /** \brief Remnant kinetic energy [MeV] */
235 /** \brief Remnant momentum, x component [MeV/c] */
237 /** \brief Remnant momentum, y component [MeV/c] */
239 /** \brief Remnant momentum, z component [MeV/c] */
241 /** \brief Remnant momentum polar angle [radians] */
243 /** \brief Remnant momentum azimuthal angle [radians] */
245 /** \brief Remnant angular momentum, x component [\f$\hbar\f$] */
247 /** \brief Remnant angular momentum, y component [\f$\hbar\f$] */
249 /** \brief Remnant angular momentum, z component [\f$\hbar\f$] */
251 /** \brief Projectile particle type */
253 /** \brief Mass number of the target nucleus */
255 /** \brief Charge number of the target nucleus */
257 /** \brief Strangeness number of the target nucleus */
259 /** \brief Mass number of the projectile nucleus */
261 /** \brief Charge number of the projectile nucleus */
263 /** \brief Strangeness number of the projectile nucleus */
265 /** \brief Projectile kinetic energy given as input */
267 /** \brief Impact parameter [fm] */
269 /** \brief Number of accepted two-body collisions */
271 /** \brief Cascade stopping time [fm/c] */
273 /** \brief Energy-conservation balance [MeV] */
275 /** \brief Longitudinal momentum-conservation balance [MeV/c] */
277 /** \brief Transverse momentum-conservation balance [MeV/c] */
279 /** \brief Number of cascade particles */
281 /** \brief True if the event is transparent */
283 /** \brief True if the event is a forced CN */
285 /** \brief True if the event is a nucleon absorption */
287 /** \brief True if the event is a pion absorption */
289 /** \brief Number of accepted Delta decays */
291 /** \brief Number of two-body collisions blocked by Pauli or CDPP */
293 /** \brief Number of decays blocked by Pauli or CDPP */
295 /** \brief Effective (Coulomb-distorted) impact parameter [fm] */
297 /** \brief Event involved deltas in the nucleus at the end of the cascade */
299 /** \brief Event involved sigmas in the nucleus at the end of the cascade */
301 /** \brief Event involved kaons in the nucleus at the end of the cascade */
303 /** \brief Event involved antikaons in the nucleus at the end of the cascade */
305 /** \brief Event involved lambdas in the nucleus at the end of the cascade */
307 /** \brief Event involved forced delta decays inside the nucleus */
309 /** \brief Event involved forced delta decays outside the nucleus */
311 /** \brief Event involved forced eta/omega decays outside the nucleus */
313 /** \brief Event involved forced strange absorption inside the nucleus */
315 /** \brief Event involved forced Sigma Zero decays outside the nucleus */
317 /** \brief Event involved forced antiKaon/Sigma absorption inside the nucleus */
319 /** \brief Number of forced Lambda emit out of the nucleus */
321 /** \brief Event involved forced Kaon emission */
323 /** \brief Event involved cluster decay */
325 /** \brief Time of the first collision [fm/c] */
327 /** \brief Cross section of the first collision (mb) */
329 /** \brief Position of the spectator on the first collision (fm) */
331 /** \brief Momentum of the spectator on the first collision (fm) */
333 /** \brief True if the first collision was elastic */
335 /** \brief Number of reflection avatars */
337 /** \brief Number of collision avatars */
339 /** \brief Number of decay avatars */
341 /** \brief Number of dynamical spectators that were merged back into the projectile remnant */
343 /** \brief Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. */
345 /** \brief Sequential number of the event in the event loop */
347 /** \brief Particle kinetic energy, in inverse kinematics [MeV] */
349 /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */
351 /** \brief Particle momentum polar angle, in inverse kinematics [radians] */
353
354 /** \brief Reset the EventInfo members */
355 void reset() {
356 nParticles = 0;
357 eventBias = (Float_t)0.0;
358 history.clear();
359 nRemnants = 0;
360 projectileType = 0;
361 At = 0;
362 Zt = 0;
363 St = 0;
364 Ap = 0;
365 Zp = 0;
366 Sp = 0;
367 Ep = (Float_t)0.0;
369 nCollisions = 0;
370 stoppingTime = (Float_t)0.0;
371 EBalance = (Float_t)0.0;
372 pLongBalance = (Float_t)0.0;
373 pTransBalance = (Float_t)0.0;
375 transparent = false;
376 forcedCompoundNucleus = false;
377 nucleonAbsorption = false;
378 pionAbsorption = false;
379 nDecays = 0;
381 nBlockedDecays = 0;
383 deltasInside = false;
384 sigmasInside = false;
385 kaonsInside = false;
386 antikaonsInside = false;
387 lambdasInside = false;
388 forcedDeltasInside = false;
389 forcedDeltasOutside = false;
392 forcedSigmaOutside = false;
393 forcedStrangeInside = false;
394 emitLambda = 0;
395 emitKaon = false;
396 clusterDecay = false;
404 nDecayAvatars = 0;
407 event = 0;
408
409 }
410
411 /// \brief Move a remnant to the particle array
412 void remnantToParticle(const G4int remnantIndex);
413
414 /// \brief Fill the variables describing the reaction in inverse kinematics
415 void fillInverseKinematics(const Double_t gamma);
416 };
417}
418
419#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
G4int Int_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.
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.
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].
Short_t Ap
Mass number of the projectile nucleus.
Int_t nReflectionAvatars
Number of reflection avatars.
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.
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.
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