Geant4  10.03.p03
 All Classes 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 // 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 
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 
58 namespace G4INCL {
59 #ifndef INCL_ROOT_USE
60  typedef G4int Int_t;
61  typedef short Short_t;
62  typedef G4float Float_t;
63  typedef G4double Double_t;
64  typedef G4bool Bool_t;
65 #endif
66 
67  struct EventInfo {
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  }
138 
141 
143  static const Short_t maxSizeRemnants = 10;
144 
146  static const Short_t maxSizeParticles = 1000;
147 
193  std::vector<std::string> history;
302 
304  void reset() {
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  }
347 
349  void remnantToParticle(const G4int remnantIndex);
350 
352  void fillInverseKinematics(const Double_t gamma);
353  };
354 }
355 
356 #endif /* G4INCLEVENTINFO_HH_HH */
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].
void fillInverseKinematics(const Double_t gamma)
Fill the variables describing the reaction in inverse kinematics.
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:89
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 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.
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 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.
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.
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.
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.