Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4INCLEventInfo.hh
이 파일의 문서화 페이지로 가기
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  St(0),
75  Ap(0),
76  Zp(0),
77  Sp(0),
78  Ep((Float_t)0.0),
79  eventBias((Float_t)0.0),
81  nCollisions(0),
82  stoppingTime((Float_t)0.0),
83  EBalance((Float_t)0.0),
84  pLongBalance((Float_t)0.0),
85  pTransBalance((Float_t)0.0),
87  transparent(false),
88  forcedCompoundNucleus(false),
89  nucleonAbsorption(false),
90  pionAbsorption(false),
91  nDecays(0),
93  nBlockedDecays(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);
134  std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0);
135  std::fill_n(ARem, maxSizeRemnants, 0);
136  std::fill_n(ZRem, maxSizeRemnants, 0);
137  std::fill_n(SRem, maxSizeRemnants, 0);
138  std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0);
139  std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0);
140  std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0);
141  std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0);
142  std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0);
143  std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0);
144  std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0);
145  std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0);
146  std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0);
147  std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0);
148  std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0);
149  std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0);
150  std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0);
151  std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0);
152  }
153 
156 
158  static const Short_t maxSizeRemnants = 10;
159 
161  static const Short_t maxSizeParticles = 1000;
162 
212  std::vector<std::string> history;
347 
349  void reset() {
350  nParticles = 0;
351  history.clear();
352  nRemnants = 0;
353  projectileType = 0;
354  At = 0;
355  Zt = 0;
356  St = 0;
357  Ap = 0;
358  Zp = 0;
359  Sp = 0;
360  Ep = (Float_t)0.0;
361  eventBias = (Float_t)0.0;
362  impactParameter = (Float_t)0.0;
363  nCollisions = 0;
364  stoppingTime = (Float_t)0.0;
365  EBalance = (Float_t)0.0;
366  pLongBalance = (Float_t)0.0;
367  pTransBalance = (Float_t)0.0;
368  nCascadeParticles = 0;
369  transparent = false;
370  forcedCompoundNucleus = false;
371  nucleonAbsorption = false;
372  pionAbsorption = false;
373  nDecays = 0;
374  nBlockedCollisions = 0;
375  nBlockedDecays = 0;
377  deltasInside = false;
378  sigmasInside = false;
379  kaonsInside = false;
380  antikaonsInside = false;
381  lambdasInside = false;
382  forcedDeltasInside = false;
383  forcedDeltasOutside = false;
385  absorbedStrangeParticle = false;
386  forcedSigmaOutside = false;
387  forcedStrangeInside = false;
388  emitLambda = 0;
389  emitKaon = false;
390  clusterDecay = false;
395  firstCollisionIsElastic = false;
396  nReflectionAvatars = 0;
397  nCollisionAvatars = 0;
398  nDecayAvatars = 0;
401  event = 0;
402 
403  }
404 
406  void remnantToParticle(const G4int remnantIndex);
407 
409  void fillInverseKinematics(const Double_t gamma);
410  };
411 }
412 
413 #endif /* G4INCLEVENTINFO_HH_HH */
Int_t event
Sequential number of the event in the event loop.
void remnantToParticle(const G4int remnantIndex)
Move a remnant to the particle array.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
G4double Double_t
Int_t nDecayAvatars
Number of decay avatars.
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
Float_t JRem[maxSizeRemnants]
Remnant spin [ ].
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
Bool_t kaonsInside
Event involved kaons in the nucleus at the end of the cascade.
Float_t firstCollisionSpectatorPosition
Position of the spectator on the first collision (fm)
Short_t St
Strangeness number of the target nucleus.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Int_t nReflectionAvatars
Number of reflection avatars.
std::vector< std::string > history
History of the particle.
void fillInverseKinematics(const Double_t gamma)
Fill the variables describing the reaction in inverse kinematics.
float G4float
Definition: G4Types.hh:77
void reset()
Reset the EventInfo members.
Short_t nParticles
Number of particles in the final state.
Short_t Z[maxSizeParticles]
Particle charge number.
Short_t nRemnants
Number of remnants.
Bool_t sigmasInside
Event involved sigmas in the nucleus at the end of the cascade.
short Short_t
Float_t Ep
Projectile kinetic energy given as input.
static const Short_t maxSizeParticles
Maximum array size for emitted particles.
static const Short_t maxSizeRemnants
Maximum array size for remnants.
#define G4ThreadLocal
Definition: tls.hh:69
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
Short_t origin[maxSizeParticles]
Origin of the particle.
Bool_t forcedSigmaOutside
Event involved forced Sigma Zero decays outside the nucleus.
Float_t pLongBalance
Longitudinal momentum-conservation balance [MeV/c].
G4float Float_t
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant. ...
Bool_t forcedCompoundNucleus
True if the event is a forced CN.
Int_t nCollisionAvatars
Number of collision avatars.
Float_t phiRem[maxSizeRemnants]
Remnant momentum azimuthal angle [radians].
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Float_t firstCollisionTime
Time of the first collision [fm/c].
Bool_t pionAbsorption
True if the event is a pion absorption.
G4int Int_t
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
Float_t impactParameter
Impact parameter [fm].
Float_t EBalance
Energy-conservation balance [MeV].
Short_t Sp
Strangeness number of the projectile nucleus.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
Float_t firstCollisionXSec
Cross section of the first collision (mb)
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t Zt
Charge number of the target nucleus.
Float_t thetaRem[maxSizeRemnants]
Remnant momentum polar angle [radians].
Bool_t forcedDeltasOutside
Event involved forced delta decays outside the nucleus.
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Bool_t antikaonsInside
Event involved antikaons in the nucleus at the end of the cascade.
Short_t S[maxSizeParticles]
Particle strangeness number.
Float_t eventBias
Event bias.
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
Bool_t forcedStrangeInside
Event involved forced antiKaon/Sigma absorbtion inside the nucleus.
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
static G4ThreadLocal Int_t eventNumber
Number of the event.
G4bool Bool_t
Bool_t emitKaon
Event involved forced Kaon emission.
int G4int
Definition: G4Types.hh:78
Bool_t forcedPionResonancesOutside
Event involved forced eta/omega decays outside the nucleus.
Int_t PDGCode[maxSizeParticles]
PDG numbering of the particles.
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
Short_t nCascadeParticles
Number of cascade particles.
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
Int_t emitLambda
Number of forced Lambda emit out of the nucleus.
Float_t EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
Int_t nCollisions
Number of accepted two-body collisions.
Bool_t deltasInside
Event involved deltas in the nucleus at the end of the cascade.
Short_t Ap
Mass number of the projectile nucleus.
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Float_t ParticleBias[maxSizeParticles]
Particle weight due to the bias.
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Short_t A[maxSizeParticles]
Particle mass number.
Float_t stoppingTime
Cascade stopping time [fm/c].
Short_t Zp
Charge number of the projectile nucleus.
Int_t projectileType
Projectile particle type.
Bool_t lambdasInside
Event involved lambdas in the nucleus at the end of the cascade.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Bool_t absorbedStrangeParticle
Event involved forced strange absorbtion inside the nucleus.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Int_t nBlockedCollisions
Number of two-body collisions blocked by Pauli or CDPP.
Short_t SRem[maxSizeRemnants]
Remnant strangeness number.
Short_t At
Mass number of the target nucleus.
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Bool_t clusterDecay
Event involved cluster decay.
Bool_t transparent
True if the event is transparent.
Int_t nDecays
Number of accepted Delta decays.
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm)