Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4INCL::EventInfo Struct Reference

#include <G4INCLEventInfo.hh>

Public Member Functions

 EventInfo ()
 
void reset ()
 Reset the EventInfo members. More...
 
void remnantToParticle (const G4int remnantIndex)
 Move a remnant to the particle array. More...
 
void fillInverseKinematics (const Double_t gamma)
 Fill the variables describing the reaction in inverse kinematics. More...
 

Public Attributes

Short_t nParticles
 Number of particles in the final state. More...
 
Short_t A [maxSizeParticles]
 Particle mass number. More...
 
Short_t Z [maxSizeParticles]
 Particle charge number. More...
 
Int_t PDGCode [maxSizeParticles]
 PDG numbering of the particles. More...
 
Float_t EKin [maxSizeParticles]
 Particle kinetic energy [MeV]. More...
 
Float_t px [maxSizeParticles]
 Particle momentum, x component [MeV/c]. More...
 
Float_t py [maxSizeParticles]
 Particle momentum, y component [MeV/c]. More...
 
Float_t pz [maxSizeParticles]
 Particle momentum, z component [MeV/c]. More...
 
Float_t theta [maxSizeParticles]
 Particle momentum polar angle [radians]. More...
 
Float_t phi [maxSizeParticles]
 Particle momentum azimuthal angle [radians]. More...
 
Short_t origin [maxSizeParticles]
 Origin of the particle. More...
 
Float_t emissionTime [maxSizeParticles]
 Emission time [fm/c]. More...
 
std::vector< std::string > history
 History of the particle. More...
 
Short_t nRemnants
 Number of remnants. More...
 
Short_t ARem [maxSizeRemnants]
 Remnant mass number. More...
 
Short_t ZRem [maxSizeRemnants]
 Remnant charge number. More...
 
Float_t EStarRem [maxSizeRemnants]
 Remnant excitation energy [MeV]. More...
 
Float_t JRem [maxSizeRemnants]
 Remnant spin [ $\hbar$]. More...
 
Float_t EKinRem [maxSizeRemnants]
 Remnant kinetic energy [MeV]. More...
 
Float_t pxRem [maxSizeRemnants]
 Remnant momentum, x component [MeV/c]. More...
 
Float_t pyRem [maxSizeRemnants]
 Remnant momentum, y component [MeV/c]. More...
 
Float_t pzRem [maxSizeRemnants]
 Remnant momentum, z component [MeV/c]. More...
 
Float_t thetaRem [maxSizeRemnants]
 Remnant momentum polar angle [radians]. More...
 
Float_t phiRem [maxSizeRemnants]
 Remnant momentum azimuthal angle [radians]. More...
 
Float_t jxRem [maxSizeRemnants]
 Remnant angular momentum, x component [ $\hbar$]. More...
 
Float_t jyRem [maxSizeRemnants]
 Remnant angular momentum, y component [ $\hbar$]. More...
 
Float_t jzRem [maxSizeRemnants]
 Remnant angular momentum, z component [ $\hbar$]. More...
 
Int_t projectileType
 Projectile particle type. More...
 
Short_t At
 Mass number of the target nucleus. More...
 
Short_t Zt
 Charge number of the target nucleus. More...
 
Short_t Ap
 Mass number of the projectile nucleus. More...
 
Short_t Zp
 Charge number of the projectile nucleus. More...
 
Float_t Ep
 Projectile kinetic energy given as input. More...
 
Float_t impactParameter
 Impact parameter [fm]. More...
 
Int_t nCollisions
 Number of accepted two-body collisions. More...
 
Float_t stoppingTime
 Cascade stopping time [fm/c]. More...
 
Float_t EBalance
 Energy-conservation balance [MeV]. More...
 
Float_t pLongBalance
 Longitudinal momentum-conservation balance [MeV/c]. More...
 
Float_t pTransBalance
 Transverse momentum-conservation balance [MeV/c]. More...
 
Short_t nCascadeParticles
 Number of cascade particles. More...
 
Bool_t transparent
 True if the event is transparent. More...
 
Bool_t forcedCompoundNucleus
 True if the event is a forced CN. More...
 
Bool_t nucleonAbsorption
 True if the event is a nucleon absorption. More...
 
Bool_t pionAbsorption
 True if the event is a pion absorption. More...
 
Int_t nDecays
 Number of accepted Delta decays. More...
 
Int_t nBlockedCollisions
 Number of two-body collisions blocked by Pauli or CDPP. More...
 
Int_t nBlockedDecays
 Number of decays blocked by Pauli or CDPP. More...
 
Float_t effectiveImpactParameter
 Effective (Coulomb-distorted) impact parameter [fm]. More...
 
Bool_t deltasInside
 Event involved deltas in the nucleus at the end of the cascade. More...
 
Bool_t forcedDeltasInside
 Event involved forced delta decays inside the nucleus. More...
 
Bool_t forcedDeltasOutside
 Event involved forced delta decays outside the nucleus. More...
 
Bool_t forcedPionResonancesOutside
 Event involved forced eta/omega decays outside the nucleus. More...
 
Bool_t clusterDecay
 Event involved cluster decay. More...
 
Float_t firstCollisionTime
 Time of the first collision [fm/c]. More...
 
Float_t firstCollisionXSec
 Cross section of the first collision (mb) More...
 
Float_t firstCollisionSpectatorPosition
 Position of the spectator on the first collision (fm) More...
 
Float_t firstCollisionSpectatorMomentum
 Momentum of the spectator on the first collision (fm) More...
 
Bool_t firstCollisionIsElastic
 True if the first collision was elastic. More...
 
Int_t nReflectionAvatars
 Number of reflection avatars. More...
 
Int_t nCollisionAvatars
 Number of collision avatars. More...
 
Int_t nDecayAvatars
 Number of decay avatars. More...
 
Int_t nUnmergedSpectators
 Number of dynamical spectators that were merged back into the projectile remnant. More...
 
Int_t nEnergyViolationInteraction
 Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. More...
 
Int_t event
 Sequential number of the event in the event loop. More...
 
Float_t EKinPrime [maxSizeParticles]
 Particle kinetic energy, in inverse kinematics [MeV]. More...
 
Float_t pzPrime [maxSizeParticles]
 Particle momentum, z component, in inverse kinematics [MeV/c]. More...
 
Float_t thetaPrime [maxSizeParticles]
 Particle momentum polar angle, in inverse kinematics [radians]. More...
 

Static Public Attributes

static G4ThreadLocal Int_t eventNumber = 0
 Number of the event. More...
 
static const Short_t maxSizeRemnants = 10
 Maximum array size for remnants. More...
 
static const Short_t maxSizeParticles = 1000
 Maximum array size for emitted particles. More...
 

Detailed Description

Definition at line 67 of file G4INCLEventInfo.hh.

Constructor & Destructor Documentation

G4INCL::EventInfo::EventInfo ( )
inline

Definition at line 68 of file G4INCLEventInfo.hh.

68  :
69  nParticles(0),
70  nRemnants(0),
71  projectileType(0),
72  At(0),
73  Zt(0),
74  Ap(0),
75  Zp(0),
76  Ep((Float_t)0.0),
78  nCollisions(0),
79  stoppingTime((Float_t)0.0),
80  EBalance((Float_t)0.0),
81  pLongBalance((Float_t)0.0),
82  pTransBalance((Float_t)0.0),
84  transparent(false),
85  forcedCompoundNucleus(false),
86  nucleonAbsorption(false),
87  pionAbsorption(false),
88  nDecays(0),
90  nBlockedDecays(0),
92  deltasInside(false),
93  forcedDeltasInside(false),
94  forcedDeltasOutside(false),
96  clusterDecay(false),
104  nDecayAvatars(0),
107  event(0)
108 
109  {
110  std::fill_n(A, maxSizeParticles, 0);
111  std::fill_n(Z, maxSizeParticles, 0);
112  std::fill_n(PDGCode, maxSizeParticles, 0);
113  std::fill_n(EKin, maxSizeParticles, (Float_t)0.0);
114  std::fill_n(px, maxSizeParticles, (Float_t)0.0);
115  std::fill_n(py, maxSizeParticles, (Float_t)0.0);
116  std::fill_n(pz, maxSizeParticles, (Float_t)0.0);
117  std::fill_n(theta, maxSizeParticles, (Float_t)0.0);
118  std::fill_n(phi, maxSizeParticles, (Float_t)0.0);
119  std::fill_n(origin, maxSizeParticles, 0);
120  std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
121  std::fill_n(ARem, maxSizeRemnants, 0);
122  std::fill_n(ZRem, maxSizeRemnants, 0);
123  std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
124  std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
125  std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
126  std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
127  std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
128  std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
129  std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
130  std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
131  std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
132  std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
133  std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
134  std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
135  std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
136  std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
137  }
Short_t Zp
Charge number of the projectile nucleus.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
Short_t nRemnants
Number of remnants.
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
Float_t pLongBalance
Longitudinal momentum-conservation balance [MeV/c].
Int_t PDGCode[maxSizeParticles]
PDG numbering of the particles.
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)
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.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
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 EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
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.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t Zt
Charge number of the target nucleus.
Float_t impactParameter
Impact parameter [fm].
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
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.
Int_t event
Sequential number of the event in the event loop.
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.
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
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.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Short_t A[maxSizeParticles]
Particle mass number.
Bool_t forcedPionResonancesOutside
Event involved forced eta/omega decays outside the nucleus.
Int_t nDecays
Number of accepted Delta decays.
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
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].
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm)
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.

Member Function Documentation

void G4INCL::EventInfo::fillInverseKinematics ( const Double_t  gamma)

Fill the variables describing the reaction in inverse kinematics.

Definition at line 56 of file G4INCLEventInfo.cc.

56  {
57  const Double_t beta = std::sqrt(1.-1./(gamma*gamma));
58  for(Int_t i=0; i<nParticles; ++i) {
59  // determine the particle mass from the kinetic energy and the momentum;
60  // this ensures consistency with the masses uses by the models
61  Double_t mass;
62  if(EKin[i]>0.) {
63  mass = std::max(
64  0.5 * (px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]-EKin[i]*EKin[i]) / EKin[i],
65  0.0);
66  } else {
67  INCL_WARN("Particle with null kinetic energy in fillInverseKinematics, cannot determine its mass:\n"
68  << " A=" << A[i] << ", Z=" << Z[i] << '\n'
69  << " EKin=" << EKin[i] << ", px=" << px[i] << ", py=" << py[i] << ", pz=" << pz[i] << '\n'
70  << " Falling back to the mass from the INCL ParticleTable" << '\n');
71  mass = ParticleTable::getRealMass(A[i], Z[i]);
72  }
73 
74  const Double_t ETot = EKin[i] + mass;
75  const Double_t ETotPrime = gamma*(ETot - beta*pz[i]);
76  EKinPrime[i] = ETotPrime - mass;
77  pzPrime[i] = -gamma*(pz[i] - beta*ETot);
78  const Double_t pPrime = std::sqrt(px[i]*px[i] + py[i]*py[i] + pzPrime[i]*pzPrime[i]);
79  const Double_t cosThetaPrime = (pPrime>0.) ? (pzPrime[i]/pPrime) : 1.;
80  if(cosThetaPrime>=1.)
81  thetaPrime[i] = 0.;
82  else if(cosThetaPrime<=-1.)
83  thetaPrime[i] = 180.;
84  else
85  thetaPrime[i] = Math::toDegrees(Math::arcCos(cosThetaPrime));
86  }
87  }
G4double arcCos(const G4double x)
Calculates arccos with some tolerance on illegal arguments.
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
#define INCL_WARN(x)
G4double toDegrees(G4double radians)
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Float_t EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
Short_t Z[maxSizeParticles]
Particle charge number.
Short_t nParticles
Number of particles in the final state.
G4int Int_t
Short_t A[maxSizeParticles]
Particle mass number.
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double Double_t
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].

Here is the call graph for this function:

void G4INCL::EventInfo::remnantToParticle ( const G4int  remnantIndex)

Move a remnant to the particle array.

Definition at line 89 of file G4INCLEventInfo.cc.

89  {
90  A[nParticles] = ARem[remnantIndex];
91  Z[nParticles] = ZRem[remnantIndex];
93 
94  px[nParticles] = pxRem[remnantIndex];
95  py[nParticles] = pyRem[remnantIndex];
96  pz[nParticles] = pzRem[remnantIndex];
97 
98  const G4double plab = std::sqrt(pxRem[remnantIndex]*pxRem[remnantIndex]
99  +pyRem[remnantIndex]*pyRem[remnantIndex]
100  +pzRem[remnantIndex]*pzRem[remnantIndex]);
101  G4double pznorm = pzRem[remnantIndex]/plab;
102  if(pznorm>1.)
103  pznorm = 1.;
104  else if(pznorm<-1.)
105  pznorm = -1.;
107  phi[nParticles] = Math::toDegrees(std::atan2(pyRem[remnantIndex],pxRem[remnantIndex]));
108 
109  EKin[nParticles] = EKinRem[remnantIndex];
110  origin[nParticles] = -1; // Origin: cascade
111  history.push_back(""); // history
112  nParticles++;
113 // assert(history.size()==(unsigned int)nParticles);
114  }
G4double arcCos(const G4double x)
Calculates arccos with some tolerance on illegal arguments.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
G4double toDegrees(G4double radians)
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t stoppingTime
Cascade stopping time [fm/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Short_t Z[maxSizeParticles]
Particle charge number.
Short_t nParticles
Number of particles in the final state.
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Short_t A[maxSizeParticles]
Particle mass number.
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
double G4double
Definition: G4Types.hh:76
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Short_t origin[maxSizeParticles]
Origin of the particle.
std::vector< std::string > history
History of the particle.

Here is the call graph for this function:

void G4INCL::EventInfo::reset ( )
inline

Reset the EventInfo members.

Definition at line 304 of file G4INCLEventInfo.hh.

304  {
305  nParticles = 0;
306  history.clear();
307  nRemnants = 0;
308  projectileType = 0;
309  At = 0;
310  Zt = 0;
311  Ap = 0;
312  Zp = 0;
313  Ep = (Float_t)0.0;
314  impactParameter = (Float_t)0.0;
315  nCollisions = 0;
316  stoppingTime = (Float_t)0.0;
317  EBalance = (Float_t)0.0;
318  pLongBalance = (Float_t)0.0;
319  pTransBalance = (Float_t)0.0;
320  nCascadeParticles = 0;
321  transparent = false;
322  forcedCompoundNucleus = false;
323  nucleonAbsorption = false;
324  pionAbsorption = false;
325  nDecays = 0;
326  nBlockedCollisions = 0;
327  nBlockedDecays = 0;
329  deltasInside = false;
330  forcedDeltasInside = false;
331  forcedDeltasOutside = false;
333  clusterDecay = false;
338  firstCollisionIsElastic = false;
339  nReflectionAvatars = 0;
340  nCollisionAvatars = 0;
341  nDecayAvatars = 0;
344  event = 0;
345 
346  }
Short_t Zp
Charge number of the projectile nucleus.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
Short_t nRemnants
Number of remnants.
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].
Int_t nCollisions
Number of accepted two-body collisions.
Float_t firstCollisionXSec
Cross section of the first collision (mb)
Bool_t forcedCompoundNucleus
True if the event is a forced CN.
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.
Int_t nReflectionAvatars
Number of reflection avatars.
Float_t stoppingTime
Cascade stopping time [fm/c].
Int_t nCollisionAvatars
Number of collision avatars.
Float_t Ep
Projectile kinetic energy given as input.
Short_t At
Mass number of the target nucleus.
Short_t Zt
Charge number of the target nucleus.
Float_t impactParameter
Impact parameter [fm].
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
Float_t EBalance
Energy-conservation balance [MeV].
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.
Int_t nDecayAvatars
Number of decay avatars.
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant. ...
Float_t firstCollisionSpectatorPosition
Position of the spectator on the first collision (fm)
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
Bool_t forcedPionResonancesOutside
Event involved forced eta/omega decays outside the nucleus.
Int_t nDecays
Number of accepted Delta decays.
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
G4float Float_t
Bool_t transparent
True if the event is transparent.
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm)
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Short_t Ap
Mass number of the projectile nucleus.
Short_t nCascadeParticles
Number of cascade particles.
Int_t projectileType
Projectile particle type.
std::vector< std::string > history
History of the particle.

Member Data Documentation

Short_t G4INCL::EventInfo::A[maxSizeParticles]

Particle mass number.

Definition at line 151 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::Ap

Mass number of the projectile nucleus.

Definition at line 229 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::ARem[maxSizeRemnants]

Remnant mass number.

Definition at line 197 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::At

Mass number of the target nucleus.

Definition at line 225 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::clusterDecay

Event involved cluster decay.

Definition at line 273 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::deltasInside

Event involved deltas in the nucleus at the end of the cascade.

Definition at line 265 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::EBalance

Energy-conservation balance [MeV].

Definition at line 241 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::effectiveImpactParameter

Effective (Coulomb-distorted) impact parameter [fm].

Definition at line 263 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::EKin[maxSizeParticles]

Particle kinetic energy [MeV].

Definition at line 157 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::EKinPrime[maxSizeParticles]

Particle kinetic energy, in inverse kinematics [MeV].

Definition at line 297 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::EKinRem[maxSizeRemnants]

Remnant kinetic energy [MeV].

Definition at line 205 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::emissionTime[maxSizeParticles]

Emission time [fm/c].

Definition at line 174 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::Ep

Projectile kinetic energy given as input.

Definition at line 233 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::EStarRem[maxSizeRemnants]

Remnant excitation energy [MeV].

Definition at line 201 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::event

Sequential number of the event in the event loop.

Definition at line 295 of file G4INCLEventInfo.hh.

G4ThreadLocal Int_t G4INCL::EventInfo::eventNumber = 0
static

Number of the event.

Definition at line 140 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::firstCollisionIsElastic

True if the first collision was elastic.

Definition at line 283 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::firstCollisionSpectatorMomentum

Momentum of the spectator on the first collision (fm)

Definition at line 281 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::firstCollisionSpectatorPosition

Position of the spectator on the first collision (fm)

Definition at line 279 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::firstCollisionTime

Time of the first collision [fm/c].

Definition at line 275 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::firstCollisionXSec

Cross section of the first collision (mb)

Definition at line 277 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::forcedCompoundNucleus

True if the event is a forced CN.

Definition at line 251 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::forcedDeltasInside

Event involved forced delta decays inside the nucleus.

Definition at line 267 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::forcedDeltasOutside

Event involved forced delta decays outside the nucleus.

Definition at line 269 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::forcedPionResonancesOutside

Event involved forced eta/omega decays outside the nucleus.

Definition at line 271 of file G4INCLEventInfo.hh.

std::vector<std::string> G4INCL::EventInfo::history

History of the particle.

Condensed information about the de-excitation chain of a particle. For cascade particles, it is just an empty string. For particles arising from the de-excitation of a cascade remnant, it is a string of characters. Each character represents one or more identical steps in the de-excitation process. The currently defined possible character values and their meanings are the following:

e: evaporation product E: evaporation residue m: multifragmentation a: light partner in asymmetric fission or IMF emission A: heavy partner in asymmetric fission or IMF emission f: light partner in fission F: heavy partner in fission s: saddle-to-scission emission n: non-statistical emission (decay)

Definition at line 193 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::impactParameter

Impact parameter [fm].

Definition at line 235 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::JRem[maxSizeRemnants]

Remnant spin [ $\hbar$].

Definition at line 203 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::jxRem[maxSizeRemnants]

Remnant angular momentum, x component [ $\hbar$].

Definition at line 217 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::jyRem[maxSizeRemnants]

Remnant angular momentum, y component [ $\hbar$].

Definition at line 219 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::jzRem[maxSizeRemnants]

Remnant angular momentum, z component [ $\hbar$].

Definition at line 221 of file G4INCLEventInfo.hh.

const Short_t G4INCL::EventInfo::maxSizeParticles = 1000
static

Maximum array size for emitted particles.

Definition at line 146 of file G4INCLEventInfo.hh.

const Short_t G4INCL::EventInfo::maxSizeRemnants = 10
static

Maximum array size for remnants.

Definition at line 143 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nBlockedCollisions

Number of two-body collisions blocked by Pauli or CDPP.

Definition at line 259 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nBlockedDecays

Number of decays blocked by Pauli or CDPP.

Definition at line 261 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::nCascadeParticles

Number of cascade particles.

Definition at line 247 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nCollisionAvatars

Number of collision avatars.

Definition at line 287 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nCollisions

Number of accepted two-body collisions.

Definition at line 237 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nDecayAvatars

Number of decay avatars.

Definition at line 289 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nDecays

Number of accepted Delta decays.

Definition at line 257 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nEnergyViolationInteraction

Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution.

Definition at line 293 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::nParticles

Number of particles in the final state.

Definition at line 149 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nReflectionAvatars

Number of reflection avatars.

Definition at line 285 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::nRemnants

Number of remnants.

Definition at line 195 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::nucleonAbsorption

True if the event is a nucleon absorption.

Definition at line 253 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::nUnmergedSpectators

Number of dynamical spectators that were merged back into the projectile remnant.

Definition at line 291 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::origin[maxSizeParticles]

Origin of the particle.

Should be -1 for cascade particles, or the number of the remnant for de-excitation particles.

Definition at line 172 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::PDGCode[maxSizeParticles]

PDG numbering of the particles.

Definition at line 155 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::phi[maxSizeParticles]

Particle momentum azimuthal angle [radians].

Definition at line 167 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::phiRem[maxSizeRemnants]

Remnant momentum azimuthal angle [radians].

Definition at line 215 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::pionAbsorption

True if the event is a pion absorption.

Definition at line 255 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pLongBalance

Longitudinal momentum-conservation balance [MeV/c].

Definition at line 243 of file G4INCLEventInfo.hh.

Int_t G4INCL::EventInfo::projectileType

Projectile particle type.

Definition at line 223 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pTransBalance

Transverse momentum-conservation balance [MeV/c].

Definition at line 245 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::px[maxSizeParticles]

Particle momentum, x component [MeV/c].

Definition at line 159 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pxRem[maxSizeRemnants]

Remnant momentum, x component [MeV/c].

Definition at line 207 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::py[maxSizeParticles]

Particle momentum, y component [MeV/c].

Definition at line 161 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pyRem[maxSizeRemnants]

Remnant momentum, y component [MeV/c].

Definition at line 209 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pz[maxSizeParticles]

Particle momentum, z component [MeV/c].

Definition at line 163 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pzPrime[maxSizeParticles]

Particle momentum, z component, in inverse kinematics [MeV/c].

Definition at line 299 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::pzRem[maxSizeRemnants]

Remnant momentum, z component [MeV/c].

Definition at line 211 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::stoppingTime

Cascade stopping time [fm/c].

Definition at line 239 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::theta[maxSizeParticles]

Particle momentum polar angle [radians].

Definition at line 165 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::thetaPrime[maxSizeParticles]

Particle momentum polar angle, in inverse kinematics [radians].

Definition at line 301 of file G4INCLEventInfo.hh.

Float_t G4INCL::EventInfo::thetaRem[maxSizeRemnants]

Remnant momentum polar angle [radians].

Definition at line 213 of file G4INCLEventInfo.hh.

Bool_t G4INCL::EventInfo::transparent

True if the event is transparent.

Definition at line 249 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::Z[maxSizeParticles]

Particle charge number.

Definition at line 153 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::Zp

Charge number of the projectile nucleus.

Definition at line 231 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::ZRem[maxSizeRemnants]

Remnant charge number.

Definition at line 199 of file G4INCLEventInfo.hh.

Short_t G4INCL::EventInfo::Zt

Charge number of the target nucleus.

Definition at line 227 of file G4INCLEventInfo.hh.


The documentation for this struct was generated from the following files: