Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::EventInfo Struct Reference

#include <G4INCLEventInfo.hh>

Public Member Functions

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

Public Attributes

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

Static Public Attributes

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

Detailed Description

Definition at line 67 of file G4INCLEventInfo.hh.

Constructor & Destructor Documentation

◆ EventInfo()

G4INCL::EventInfo::EventInfo ( )
inline

Definition at line 68 of file G4INCLEventInfo.hh.

68 :
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 }
G4float Float_t
Bool_t pionAbsorption
True if the event is a pion absorption.
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].
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].
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.
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.
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.

Member Function Documentation

◆ fillInverseKinematics()

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

Fill the variables describing the reaction in inverse kinematics.

Definition at line 57 of file G4INCLEventInfo.cc.

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

◆ remnantToParticle()

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

Move a remnant to the particle array.

Definition at line 90 of file G4INCLEventInfo.cc.

90 {
91
92 INCL_DEBUG("remnantToParticle function used\n");
93
94 A[nParticles] = ARem[remnantIndex];
95 Z[nParticles] = ZRem[remnantIndex];
96 S[nParticles] = SRem[remnantIndex];
97
98 ParticleSpecies pt(A[nParticles],Z[nParticles],S[nParticles]);
99 PDGCode[nParticles] = pt.getPDGCode();
100
103
104 px[nParticles] = pxRem[remnantIndex];
105 py[nParticles] = pyRem[remnantIndex];
106 pz[nParticles] = pzRem[remnantIndex];
107
108 const G4double plab = std::sqrt(pxRem[remnantIndex]*pxRem[remnantIndex]
109 +pyRem[remnantIndex]*pyRem[remnantIndex]
110 +pzRem[remnantIndex]*pzRem[remnantIndex]);
111 G4double pznorm = pzRem[remnantIndex]/plab;
112 if(pznorm>1.)
113 pznorm = 1.;
114 else if(pznorm<-1.)
115 pznorm = -1.;
117 phi[nParticles] = Math::toDegrees(std::atan2(pyRem[remnantIndex],pxRem[remnantIndex]));
118
119 EKin[nParticles] = EKinRem[remnantIndex];
120 origin[nParticles] = -1; // Origin: cascade
121 parentResonancePDGCode[nParticles] = 0; // No parent resonance
122 parentResonanceID[nParticles] = 0; // No parent resonance
123 history.push_back(""); // history
124 nParticles++;
125// assert(history.size()==(unsigned int)nParticles);
126 }
#define INCL_DEBUG(x)
double G4double
Definition: G4Types.hh:83
static G4double getTotalBias()
General bias vector function.
std::vector< std::string > history
History of the particle.

◆ reset()

void G4INCL::EventInfo::reset ( )
inline

Reset the EventInfo members.

Definition at line 355 of file G4INCLEventInfo.hh.

355 {
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 }

Member Data Documentation

◆ A

◆ absorbedStrangeParticle

Bool_t G4INCL::EventInfo::absorbedStrangeParticle

Event involved forced strange absorption inside the nucleus.

Definition at line 314 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ antikaonsInside

Bool_t G4INCL::EventInfo::antikaonsInside

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

Definition at line 304 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ Ap

Short_t G4INCL::EventInfo::Ap

Mass number of the projectile nucleus.

Definition at line 260 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().

◆ ARem

Short_t G4INCL::EventInfo::ARem[maxSizeRemnants]

Remnant mass number.

Definition at line 224 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ At

Short_t G4INCL::EventInfo::At

Mass number of the target nucleus.

Definition at line 254 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().

◆ clusterDecay

Bool_t G4INCL::EventInfo::clusterDecay

Event involved cluster decay.

Definition at line 324 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ deltasInside

Bool_t G4INCL::EventInfo::deltasInside

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

Definition at line 298 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ EBalance

Float_t G4INCL::EventInfo::EBalance

Energy-conservation balance [MeV].

Definition at line 274 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ effectiveImpactParameter

Float_t G4INCL::EventInfo::effectiveImpactParameter

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

Definition at line 296 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ EKin

◆ EKinPrime

Float_t G4INCL::EventInfo::EKinPrime[maxSizeParticles]

Particle kinetic energy, in inverse kinematics [MeV].

Definition at line 348 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and fillInverseKinematics().

◆ EKinRem

Float_t G4INCL::EventInfo::EKinRem[maxSizeRemnants]

Remnant kinetic energy [MeV].

Definition at line 234 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ emissionTime

Float_t G4INCL::EventInfo::emissionTime[maxSizeParticles]

Emission time [fm/c].

Definition at line 201 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ emitKaon

Bool_t G4INCL::EventInfo::emitKaon

Event involved forced Kaon emission.

Definition at line 322 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ emitLambda

Int_t G4INCL::EventInfo::emitLambda

Number of forced Lambda emit out of the nucleus.

Definition at line 320 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ Ep

Float_t G4INCL::EventInfo::Ep

Projectile kinetic energy given as input.

Definition at line 266 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ EStarRem

Float_t G4INCL::EventInfo::EStarRem[maxSizeRemnants]

Remnant excitation energy [MeV].

Definition at line 230 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ event

Int_t G4INCL::EventInfo::event

Sequential number of the event in the event loop.

Definition at line 346 of file G4INCLEventInfo.hh.

◆ eventBias

Float_t G4INCL::EventInfo::eventBias

Event bias.

Definition at line 178 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ eventNumber

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

Number of the event.

Definition at line 157 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo().

◆ firstCollisionIsElastic

Bool_t G4INCL::EventInfo::firstCollisionIsElastic

True if the first collision was elastic.

Definition at line 334 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ firstCollisionSpectatorMomentum

Float_t G4INCL::EventInfo::firstCollisionSpectatorMomentum

Momentum of the spectator on the first collision (fm)

Definition at line 332 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ firstCollisionSpectatorPosition

Float_t G4INCL::EventInfo::firstCollisionSpectatorPosition

Position of the spectator on the first collision (fm)

Definition at line 330 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ firstCollisionTime

Float_t G4INCL::EventInfo::firstCollisionTime

Time of the first collision [fm/c].

Definition at line 326 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ firstCollisionXSec

Float_t G4INCL::EventInfo::firstCollisionXSec

Cross section of the first collision (mb)

Definition at line 328 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ forcedCompoundNucleus

Bool_t G4INCL::EventInfo::forcedCompoundNucleus

True if the event is a forced CN.

Definition at line 284 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ forcedDeltasInside

Bool_t G4INCL::EventInfo::forcedDeltasInside

Event involved forced delta decays inside the nucleus.

Definition at line 308 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ forcedDeltasOutside

Bool_t G4INCL::EventInfo::forcedDeltasOutside

Event involved forced delta decays outside the nucleus.

Definition at line 310 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ forcedPionResonancesOutside

Bool_t G4INCL::EventInfo::forcedPionResonancesOutside

Event involved forced eta/omega decays outside the nucleus.

Definition at line 312 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ forcedSigmaOutside

Bool_t G4INCL::EventInfo::forcedSigmaOutside

Event involved forced Sigma Zero decays outside the nucleus.

Definition at line 316 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ forcedStrangeInside

Bool_t G4INCL::EventInfo::forcedStrangeInside

Event involved forced antiKaon/Sigma absorption inside the nucleus.

Definition at line 318 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ history

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 220 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), remnantToParticle(), and reset().

◆ impactParameter

Float_t G4INCL::EventInfo::impactParameter

Impact parameter [fm].

Definition at line 268 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ JRem

Float_t G4INCL::EventInfo::JRem[maxSizeRemnants]

Remnant spin [ $\hbar$].

Definition at line 232 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ jxRem

Float_t G4INCL::EventInfo::jxRem[maxSizeRemnants]

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

Definition at line 246 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ jyRem

Float_t G4INCL::EventInfo::jyRem[maxSizeRemnants]

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

Definition at line 248 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ jzRem

Float_t G4INCL::EventInfo::jzRem[maxSizeRemnants]

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

Definition at line 250 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ kaonsInside

Bool_t G4INCL::EventInfo::kaonsInside

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

Definition at line 302 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ lambdasInside

Bool_t G4INCL::EventInfo::lambdasInside

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

Definition at line 306 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ maxSizeParticles

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

Maximum array size for emitted particles.

Definition at line 163 of file G4INCLEventInfo.hh.

Referenced by EventInfo().

◆ maxSizeRemnants

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

Maximum array size for remnants.

Definition at line 160 of file G4INCLEventInfo.hh.

Referenced by EventInfo().

◆ nBlockedCollisions

Int_t G4INCL::EventInfo::nBlockedCollisions

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

Definition at line 292 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nBlockedDecays

Int_t G4INCL::EventInfo::nBlockedDecays

Number of decays blocked by Pauli or CDPP.

Definition at line 294 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nCascadeParticles

Short_t G4INCL::EventInfo::nCascadeParticles

Number of cascade particles.

Definition at line 280 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nCollisionAvatars

Int_t G4INCL::EventInfo::nCollisionAvatars

Number of collision avatars.

Definition at line 338 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nCollisions

Int_t G4INCL::EventInfo::nCollisions

Number of accepted two-body collisions.

Definition at line 270 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nDecayAvatars

Int_t G4INCL::EventInfo::nDecayAvatars

Number of decay avatars.

Definition at line 340 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nDecays

Int_t G4INCL::EventInfo::nDecays

Number of accepted Delta decays.

Definition at line 290 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nEnergyViolationInteraction

Int_t G4INCL::EventInfo::nEnergyViolationInteraction

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

Definition at line 344 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nParticles

Short_t G4INCL::EventInfo::nParticles

Number of particles in the final state.

Definition at line 166 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), G4INCL::Nucleus::fillEventInfo(), fillInverseKinematics(), remnantToParticle(), and reset().

◆ nReflectionAvatars

Int_t G4INCL::EventInfo::nReflectionAvatars

Number of reflection avatars.

Definition at line 336 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nRemnants

Short_t G4INCL::EventInfo::nRemnants

Number of remnants.

Definition at line 222 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nucleonAbsorption

Bool_t G4INCL::EventInfo::nucleonAbsorption

True if the event is a nucleon absorption.

Definition at line 286 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ nUnmergedSpectators

Int_t G4INCL::EventInfo::nUnmergedSpectators

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

Definition at line 342 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ origin

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 195 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ parentResonanceID

Int_t G4INCL::EventInfo::parentResonanceID[maxSizeParticles]

Particle's parent resonance unique ID identifier.

Definition at line 199 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ parentResonancePDGCode

Int_t G4INCL::EventInfo::parentResonancePDGCode[maxSizeParticles]

Particle's parent resonance PDG code.

Definition at line 197 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ ParticleBias

Float_t G4INCL::EventInfo::ParticleBias[maxSizeParticles]

Particle weight due to the bias.

Definition at line 176 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ PDGCode

Int_t G4INCL::EventInfo::PDGCode[maxSizeParticles]

PDG numbering of the particles.

Definition at line 174 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ phi

Float_t G4INCL::EventInfo::phi[maxSizeParticles]

Particle momentum azimuthal angle [radians].

Definition at line 190 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ phiRem

Float_t G4INCL::EventInfo::phiRem[maxSizeRemnants]

Remnant momentum azimuthal angle [radians].

Definition at line 244 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ pionAbsorption

Bool_t G4INCL::EventInfo::pionAbsorption

True if the event is a pion absorption.

Definition at line 288 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ pLongBalance

Float_t G4INCL::EventInfo::pLongBalance

Longitudinal momentum-conservation balance [MeV/c].

Definition at line 276 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ projectileType

Int_t G4INCL::EventInfo::projectileType

Projectile particle type.

Definition at line 252 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::fillEventInfo(), and reset().

◆ pTransBalance

Float_t G4INCL::EventInfo::pTransBalance

Transverse momentum-conservation balance [MeV/c].

Definition at line 278 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ px

Float_t G4INCL::EventInfo::px[maxSizeParticles]

◆ pxRem

Float_t G4INCL::EventInfo::pxRem[maxSizeRemnants]

Remnant momentum, x component [MeV/c].

Definition at line 236 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ py

Float_t G4INCL::EventInfo::py[maxSizeParticles]

◆ pyRem

Float_t G4INCL::EventInfo::pyRem[maxSizeRemnants]

Remnant momentum, y component [MeV/c].

Definition at line 238 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ pz

Float_t G4INCL::EventInfo::pz[maxSizeParticles]

◆ pzPrime

Float_t G4INCL::EventInfo::pzPrime[maxSizeParticles]

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

Definition at line 350 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and fillInverseKinematics().

◆ pzRem

Float_t G4INCL::EventInfo::pzRem[maxSizeRemnants]

Remnant momentum, z component [MeV/c].

Definition at line 240 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ S

◆ sigmasInside

Bool_t G4INCL::EventInfo::sigmasInside

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

Definition at line 300 of file G4INCLEventInfo.hh.

Referenced by reset().

◆ Sp

Short_t G4INCL::EventInfo::Sp

Strangeness number of the projectile nucleus.

Definition at line 264 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().

◆ SRem

Short_t G4INCL::EventInfo::SRem[maxSizeRemnants]

Remnant strangeness number.

Definition at line 228 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ St

Short_t G4INCL::EventInfo::St

Strangeness number of the target nucleus.

Definition at line 258 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().

◆ stoppingTime

Float_t G4INCL::EventInfo::stoppingTime

Cascade stopping time [fm/c].

Definition at line 272 of file G4INCLEventInfo.hh.

Referenced by remnantToParticle(), and reset().

◆ theta

Float_t G4INCL::EventInfo::theta[maxSizeParticles]

Particle momentum polar angle [radians].

Definition at line 188 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ thetaPrime

Float_t G4INCL::EventInfo::thetaPrime[maxSizeParticles]

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

Definition at line 352 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and fillInverseKinematics().

◆ thetaRem

Float_t G4INCL::EventInfo::thetaRem[maxSizeRemnants]

Remnant momentum polar angle [radians].

Definition at line 242 of file G4INCLEventInfo.hh.

Referenced by EventInfo(), and G4INCL::Nucleus::fillEventInfo().

◆ transparent

Bool_t G4INCL::EventInfo::transparent

True if the event is transparent.

Definition at line 282 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), G4INCL::INCL::processEvent(), and reset().

◆ Z

◆ Zp

Short_t G4INCL::EventInfo::Zp

Charge number of the projectile nucleus.

Definition at line 262 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().

◆ ZRem

Short_t G4INCL::EventInfo::ZRem[maxSizeRemnants]

Remnant charge number.

Definition at line 226 of file G4INCLEventInfo.hh.

Referenced by G4INCLXXInterface::ApplyYourself(), EventInfo(), G4INCL::Nucleus::fillEventInfo(), and remnantToParticle().

◆ Zt

Short_t G4INCL::EventInfo::Zt

Charge number of the target nucleus.

Definition at line 256 of file G4INCLEventInfo.hh.

Referenced by G4INCL::Nucleus::getConservationBalance(), and reset().


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