Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // 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 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
37 /** \file G4INCLEventInfo.hh
38  * \brief Simple container for output of event results.
39  *
40  * Contains the results of an INCL cascade.
41  *
42  * \date 21 January 2011
43  * \author Davide Mancusi
44  */
45 
46 #ifndef G4INCLEVENTINFO_HH_HH
47 #define G4INCLEVENTINFO_HH_HH 1
48 
49 #include "G4INCLParticleType.hh"
50 #ifdef INCL_ROOT_USE
51 #include <Rtypes.h>
52 #endif
53 #include <string>
54 #include <vector>
55 #include <algorithm>
56 
57 namespace G4INCL {
58 #ifndef INCL_ROOT_USE
59  typedef G4int Int_t;
60  typedef short Short_t;
61  typedef G4float Float_t;
62  typedef G4double Double_t;
63  typedef G4bool Bool_t;
64 #endif
65 
66  struct EventInfo {
68  nParticles(0),
69  nRemnants(0),
70  projectileType(0),
71  At(0),
72  Zt(0),
73  Ap(0),
74  Zp(0),
75  Ep((Float_t)0.0),
77  nCollisions(0),
78  stoppingTime((Float_t)0.0),
79  EBalance((Float_t)0.0),
80  pLongBalance((Float_t)0.0),
81  pTransBalance((Float_t)0.0),
83  transparent(false),
84  forcedCompoundNucleus(false),
85  nucleonAbsorption(false),
86  pionAbsorption(false),
87  nDecays(0),
89  nBlockedDecays(0),
91  deltasInside(false),
92  forcedDeltasInside(false),
93  forcedDeltasOutside(false),
94  clusterDecay(false),
102  nDecayAvatars(0),
105 #ifdef INCL_INVERSE_KINEMATICS
106 
107 #endif
108  {
109  std::fill_n(A, maxSizeParticles, 0);
110  std::fill_n(Z, maxSizeParticles, 0);
111  std::fill_n(EKin, maxSizeParticles, (Float_t)0.0);
112  std::fill_n(px, maxSizeParticles, (Float_t)0.0);
113  std::fill_n(py, maxSizeParticles, (Float_t)0.0);
114  std::fill_n(pz, maxSizeParticles, (Float_t)0.0);
115  std::fill_n(theta, maxSizeParticles, (Float_t)0.0);
116  std::fill_n(phi, maxSizeParticles, (Float_t)0.0);
117  std::fill_n(origin, maxSizeParticles, 0);
118  std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
119  std::fill_n(ARem, maxSizeRemnants, 0);
120  std::fill_n(ZRem, maxSizeRemnants, 0);
121  std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
122  std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
123  std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
124  std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
125  std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
126  std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
127  std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
128  std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
129  std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
130  std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
131  std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
132 #ifdef INCL_INVERSE_KINEMATICS
133  std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
134  std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
135  std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
136 #endif
137  }
138 
139  /** \brief Number of the event */
141 
142  /** \brief Maximum array size for remnants */
143  static const Short_t maxSizeRemnants = 10;
144 
145  /** \brief Maximum array size for emitted particles */
146  static const Short_t maxSizeParticles = 1000;
147 
148  /** \brief Number of particles in the final state */
150  /** \brief Particle mass number */
152  /** \brief Particle charge number */
154  /** \brief Particle kinetic energy [MeV] */
156  /** \brief Particle momentum, x component [MeV/c] */
158  /** \brief Particle momentum, y component [MeV/c] */
160  /** \brief Particle momentum, z component [MeV/c] */
162  /** \brief Particle momentum polar angle [radians] */
164  /** \brief Particle momentum azimuthal angle [radians] */
166  /** \brief Origin of the particle
167  *
168  * Should be -1 for cascade particles, or the number of the remnant for
169  * de-excitation particles. */
171  /** \brief Emission time [fm/c] */
173  /** \brief History of the particle
174  *
175  * Condensed information about the de-excitation chain of a particle. For
176  * cascade particles, it is just an empty string. For particles arising
177  * from the de-excitation of a cascade remnant, it is a string of
178  * characters. Each character represents one or more identical steps in
179  * the de-excitation process. The currently defined possible character
180  * values and their meanings are the following:
181  *
182  * e: evaporation product
183  * E: evaporation residue
184  * m: multifragmentation
185  * a: light partner in asymmetric fission or IMF emission
186  * A: heavy partner in asymmetric fission or IMF emission
187  * f: light partner in fission
188  * F: heavy partner in fission
189  * s: saddle-to-scission emission
190  * n: non-statistical emission (decay) */
191  std::vector<std::string> history;
192  /** \brief Number of remnants */
194  /** \brief Remnant mass number */
196  /** \brief Remnant charge number */
198  /** \brief Remnant excitation energy [MeV] */
200  /** \brief Remnant spin [\f$\hbar\f$] */
202  /** \brief Remnant kinetic energy [MeV] */
204  /** \brief Remnant momentum, x component [MeV/c] */
206  /** \brief Remnant momentum, y component [MeV/c] */
208  /** \brief Remnant momentum, z component [MeV/c] */
210  /** \brief Remnant momentum polar angle [radians] */
212  /** \brief Remnant momentum azimuthal angle [radians] */
214  /** \brief Remnant angular momentum, x component [\f$\hbar\f$] */
216  /** \brief Remnant angular momentum, y component [\f$\hbar\f$] */
218  /** \brief Remnant angular momentum, z component [\f$\hbar\f$] */
220  /** \brief Projectile particle type */
222  /** \brief Mass number of the target nucleus */
224  /** \brief Charge number of the target nucleus */
226  /** \brief Mass number of the projectile nucleus */
228  /** \brief Charge number of the projectile nucleus */
230  /** \brief Projectile kinetic energy given as input */
232  /** \brief Impact parameter [fm] */
234  /** \brief Number of accepted two-body collisions */
236  /** \brief Cascade stopping time [fm/c] */
238  /** \brief Energy-conservation balance [MeV] */
240  /** \brief Longitudinal momentum-conservation balance [MeV/c] */
242  /** \brief Transverse momentum-conservation balance [MeV/c] */
244  /** \brief Number of cascade particles */
246  /** \brief True if the event is transparent */
248  /** \brief True if the event is a forced CN */
250  /** \brief True if the event is a nucleon absorption */
252  /** \brief True if the event is a pion absorption */
254  /** \brief Number of accepted Delta decays */
256  /** \brief Number of two-body collisions blocked by Pauli or CDPP */
258  /** \brief Number of decays blocked by Pauli or CDPP */
260  /** \brief Effective (Coulomb-distorted) impact parameter [fm] */
262  /** \brief Event involved deltas in the nucleus at the end of the cascade */
264  /** \brief Event involved forced delta decays inside the nucleus */
266  /** \brief Event involved forced delta decays outside the nucleus */
268  /** \brief Event involved cluster decay */
270  /** \brief Time of the first collision [fm/c] */
272  /** \brief Cross section of the first collision (mb) */
274  /** \brief Position of the spectator on the first collision (fm) */
276  /** \brief Momentum of the spectator on the first collision (fm) */
278  /** \brief True if the first collision was elastic */
280  /** \brief Number of reflection avatars */
282  /** \brief Number of collision avatars */
284  /** \brief Number of decay avatars */
286  /** \brief Number of dynamical spectators that were merged back into the projectile remnant */
288  /** \brief Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. */
290 
291 #ifdef INCL_INVERSE_KINEMATICS
292  /** \brief Particle kinetic energy, in inverse kinematics [MeV] */
293  Float_t EKinPrime[maxSizeParticles];
294  /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */
295  Float_t pzPrime[maxSizeParticles];
296  /** \brief Particle momentum polar angle, in inverse kinematics [radians] */
297  Float_t thetaPrime[maxSizeParticles];
298 #endif // INCL_INVERSE_KINEMATICS
299 
300  /** \brief Reset the EventInfo members */
301  void reset() {
302  nParticles = 0;
303  history.clear();
304  nRemnants = 0;
305  projectileType = 0;
306  At = 0;
307  Zt = 0;
308  Ap = 0;
309  Zp = 0;
310  Ep = (Float_t)0.0;
311  impactParameter = (Float_t)0.0;
312  nCollisions = 0;
313  stoppingTime = (Float_t)0.0;
314  EBalance = (Float_t)0.0;
315  pLongBalance = (Float_t)0.0;
316  pTransBalance = (Float_t)0.0;
317  nCascadeParticles = 0;
318  transparent = false;
319  forcedCompoundNucleus = false;
320  nucleonAbsorption = false;
321  pionAbsorption = false;
322  nDecays = 0;
323  nBlockedCollisions = 0;
324  nBlockedDecays = 0;
326  deltasInside = false;
327  forcedDeltasInside = false;
328  forcedDeltasOutside = false;
329  clusterDecay = false;
334  firstCollisionIsElastic = false;
335  nReflectionAvatars = 0;
336  nCollisionAvatars = 0;
337  nDecayAvatars = 0;
340 #ifdef INCL_INVERSE_KINEMATICS
341 
342 #endif
343  }
344 
345  /// \brief Move a remnant to the particle array
346  void remnantToParticle(const G4int remnantIndex);
347 
348 #ifdef INCL_INVERSE_KINEMATICS
349  /// \brief Fill the variables describing the reaction in inverse kinematics
350  void fillInverseKinematics(const Double_t gamma);
351 #endif // INCL_INVERSE_KINEMATICS
352  };
353 }
354 
355 #endif /* G4INCLEVENTINFO_HH_HH */
Short_t Zp
Charge number of the projectile nucleus.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
Float_t pLongBalance
Longitudinal momentum-conservation balance [MeV/c].
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Float_t firstCollisionTime
Time of the first collision [fm/c].
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Int_t nCollisions
Number of accepted two-body collisions.
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
Float_t firstCollisionXSec
Cross section of the first collision (mb)
float G4float
Definition: G4Types.hh:77
Bool_t forcedCompoundNucleus
True if the event is a forced CN.
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
Bool_t forcedDeltasOutside
Event involved forced delta decays outside the nucleus.
Int_t nBlockedCollisions
Number of two-body collisions blocked by Pauli or CDPP.
Bool_t pionAbsorption
True if the event is a pion absorption.
Bool_t clusterDecay
Event involved cluster decay.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Float_t JRem[maxSizeRemnants]
Remnant spin [ ].
Int_t nReflectionAvatars
Number of reflection avatars.
#define G4ThreadLocal
Definition: tls.hh:52
short Short_t
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
int G4int
Definition: G4Types.hh:78
Float_t stoppingTime
Cascade stopping time [fm/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Int_t nCollisionAvatars
Number of collision avatars.
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Float_t Ep
Projectile kinetic energy given as input.
Short_t At
Mass number of the target nucleus.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
static G4ThreadLocal Int_t eventNumber
Number of the event.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t Zt
Charge number of the target nucleus.
Float_t impactParameter
Impact parameter [fm].
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
Float_t EBalance
Energy-conservation balance [MeV].
Short_t Z[maxSizeParticles]
Particle charge number.
Short_t nParticles
Number of particles in the final state.
Bool_t deltasInside
Event involved deltas in the nucleus at the end of the cascade.
G4int Int_t
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
bool G4bool
Definition: G4Types.hh:79
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
Int_t nDecayAvatars
Number of decay avatars.
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant. ...
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
Float_t firstCollisionSpectatorPosition
Position of the spectator on the first collision (fm)
Float_t phiRem[maxSizeRemnants]
Remnant momentum azimuthal angle [radians].
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
void remnantToParticle(const G4int remnantIndex)
Move a remnant to the particle array.
Int_t nRemnants
Number of remnants.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Short_t A[maxSizeParticles]
Particle mass number.
Int_t nDecays
Number of accepted Delta decays.
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
void reset()
Reset the EventInfo members.
Float_t thetaRem[maxSizeRemnants]
Remnant momentum polar angle [radians].
G4float Float_t
static const Short_t maxSizeRemnants
Maximum array size for remnants.
static const Short_t maxSizeParticles
Maximum array size for emitted particles.
Bool_t transparent
True if the event is transparent.
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
G4double Double_t
double G4double
Definition: G4Types.hh:76
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm)
G4bool Bool_t
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Short_t Ap
Mass number of the projectile nucleus.
Short_t nCascadeParticles
Number of cascade particles.
Short_t origin[maxSizeParticles]
Origin of the particle.
Int_t projectileType
Projectile particle type.
std::vector< std::string > history
History of the particle.