Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Nucleus.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// original by H.P. Wellisch
27// modified by J.L. Chuma, TRIUMF, 19-Nov-1996
28// last modified: 27-Mar-1997
29// Chr. Volcker, 10-Nov-1997: new methods and class variables.
30// M.G. Pia, 2 Oct 1998: modified GetFermiMomentum (original design was
31// the source of memory leaks)
32// G.Folger, spring 2010: add integer A/Z interface
33// A. Ribon, autumn 2021: extended to hypernuclei
34
35#ifndef G4Nucleus_h
36#define G4Nucleus_h 1
37// Class Description
38// This class knows how to describe a nucleus;
39// to be used in your physics implementation (not physics list) in case you need this physics.
40// Class Description - End
41
42
43#include "globals.hh"
44#include "G4ThreeVector.hh"
45#include "G4ParticleTypes.hh"
46#include "G4ReactionProduct.hh"
47#include "G4DynamicParticle.hh"
49#include "Randomize.hh"
50
52{
53 public:
54
55 G4Nucleus();
56 G4Nucleus(const G4double A, const G4double Z, const G4int numberOfLambdas = 0);
57 G4Nucleus(const G4int A, const G4int Z, const G4int numberOfLambdas = 0);
58 G4Nucleus(const G4Material* aMaterial);
59
60 ~G4Nucleus();
61
62 G4Nucleus(const G4Nucleus&) = default;
63 G4Nucleus(G4Nucleus&&) = default;
64 G4Nucleus& operator = (const G4Nucleus&) = default;
66
67 inline G4bool operator==( const G4Nucleus &right ) const
68 { return ( this == (G4Nucleus *) &right ); }
69
70 inline G4bool operator!=( const G4Nucleus &right ) const
71 { return ( this != (G4Nucleus *) &right ); }
72
73 void ChooseParameters( const G4Material *aMaterial );
74
75 void SetParameters( const G4double A, const G4double Z, const G4int numberOfLambdas = 0 );
76 void SetParameters( const G4int A, const G4int Z, const G4int numberOfLambdas = 0 );
77
78 inline G4int GetA_asInt() const
79 { return theA; }
80
81 inline G4int GetN_asInt() const
82 { return theA-theZ-theL; }
83
84 inline G4int GetZ_asInt() const
85 { return theZ; }
86
87 inline G4int GetL() const // Number of Lambdas (in the case of a hypernucleus)
88 { return theL; }
89
90 inline const G4Isotope* GetIsotope()
91 { return fIsotope; }
92
93 inline void SetIsotope(const G4Isotope* iso)
94 {
95 fIsotope = iso;
96 if(iso) {
97 theZ = iso->GetZ();
98 theA = iso->GetN();
99 theL = 0;
100 aEff = theA;
101 zEff = theZ;
102 }
103 }
104
106
107 G4double AtomicMass( const G4double A, const G4double Z, const G4int numberOfLambdas = 0 ) const;
108 G4double AtomicMass( const G4int A, const G4int Z, const G4int numberOfLambdas = 0 ) const;
109
110 G4double GetThermalPz( const G4double mass, const G4double temp ) const;
111
113
115
116 void DoKinematicsOfThermalNucleus(const G4double mu, const G4double vT_norm, const G4ThreeVector& aVelocity,
117 G4ReactionProduct& result) const;
118
119 G4double Cinema( G4double kineticEnergy );
120
121 G4double EvaporationEffects( G4double kineticEnergy );
122
124
126 { return pnBlackTrackEnergy; }
127
129 { return dtaBlackTrackEnergy; }
130
132 { return pnBlackTrackEnergyfromAnnihilation; }
133
135 { return dtaBlackTrackEnergyfromAnnihilation; }
136
137// ****************** methods introduced by ChV ***********************
138 // return fermi momentum
140
141/*
142 // return particle to be absorbed.
143 G4DynamicParticle* ReturnAbsorbingParticle(G4double weight);
144*/
145
146 // final nucleus fragmentation. Return List of particles
147 // which should be used for further tracking.
149
150
151 // excitation Energy...
152 void AddExcitationEnergy(G4double anEnergy);
153
154
155 // momentum of absorbed Particles ..
156 void AddMomentum(const G4ThreeVector aMomentum);
157
158 // return excitation Energy
159 G4double GetEnergyDeposit() {return excitationEnergy; }
160
161
162
163// ****************************** end ChV ******************************
164
165
166 private:
167
168 G4int theA;
169 G4int theZ;
170 G4int theL; // Number of Lambdas (in the case of hypernucleus)
171 G4double aEff; // effective atomic weight
172 G4double zEff; // effective atomic number
173
174 const G4Isotope* fIsotope;
175
176 G4double pnBlackTrackEnergy; // the kinetic energy available for
177 // proton/neutron black track particles
178 G4double dtaBlackTrackEnergy; // the kinetic energy available for
179 // deuteron/triton/alpha particles
180 G4double pnBlackTrackEnergyfromAnnihilation;
181 // kinetic energy available for proton/neutron black
182 // track particles based on baryon annihilation
183 G4double dtaBlackTrackEnergyfromAnnihilation;
184 // kinetic energy available for deuteron/triton/alpha
185 // black track particles based on baryon annihilation
186
187
188// ************************** member variables by ChV *******************
189 // Excitation Energy leading to evaporation or deexcitation.
190 G4double excitationEnergy;
191
192 // Momentum, accumulated by absorbing Particles
193 G4ThreeVector momentum;
194
195 // Fermi Gas model: at present, we assume constant nucleon density for all
196 // nuclei. The radius of a nucleon is taken to be 1 fm.
197 // see for example S.Fl"ugge, Encyclopedia of Physics, Vol XXXIX,
198 // Structure of Atomic Nuclei (Berlin-Gottingen-Heidelberg, 1957) page 426.
199
200 // maximum momentum possible from fermi gas model:
201 G4double fermiMomentum;
202 G4double theTemp; // temperature
203// ****************************** end ChV ******************************
204
205 };
206
207#endif
208
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4int GetZ() const
Definition G4Isotope.hh:80
G4int GetN() const
Definition G4Isotope.hh:83
void AddExcitationEnergy(G4double anEnergy)
Definition G4Nucleus.cc:578
G4double GetThermalPz(const G4double mass, const G4double temp) const
Definition G4Nucleus.cc:393
G4int GetA_asInt() const
Definition G4Nucleus.hh:78
G4int GetZ_asInt() const
Definition G4Nucleus.hh:84
G4double EvaporationEffects(G4double kineticEnergy)
Definition G4Nucleus.cc:404
void ChooseParameters(const G4Material *aMaterial)
Definition G4Nucleus.cc:277
G4double GetAnnihilationPNBlackTrackEnergy() const
Definition G4Nucleus.hh:131
G4int GetL() const
Definition G4Nucleus.hh:87
G4double AtomicMass(const G4double A, const G4double Z, const G4int numberOfLambdas=0) const
Definition G4Nucleus.cc:369
G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg)
Definition G4Nucleus.cc:465
const G4Isotope * GetIsotope()
Definition G4Nucleus.hh:90
G4double GetPNBlackTrackEnergy() const
Definition G4Nucleus.hh:125
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition G4Nucleus.cc:319
G4double Cinema(G4double kineticEnergy)
Definition G4Nucleus.cc:511
G4double GetAnnihilationDTABlackTrackEnergy() const
Definition G4Nucleus.hh:134
G4bool operator!=(const G4Nucleus &right) const
Definition G4Nucleus.hh:70
G4int GetN_asInt() const
Definition G4Nucleus.hh:81
G4DynamicParticle * ReturnTargetParticle() const
Definition G4Nucleus.cc:352
G4double GetEnergyDeposit()
Definition G4Nucleus.hh:159
G4ReactionProductVector * Fragmentate()
Definition G4Nucleus.cc:565
G4Nucleus(const G4Nucleus &)=default
void SetIsotope(const G4Isotope *iso)
Definition G4Nucleus.hh:93
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition G4Nucleus.cc:119
G4double GetDTABlackTrackEnergy() const
Definition G4Nucleus.hh:128
void AddMomentum(const G4ThreeVector aMomentum)
Definition G4Nucleus.cc:572
void DoKinematicsOfThermalNucleus(const G4double mu, const G4double vT_norm, const G4ThreeVector &aVelocity, G4ReactionProduct &result) const
Definition G4Nucleus.cc:190
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition G4Nucleus.cc:248
G4Nucleus(G4Nucleus &&)=default
G4Nucleus & operator=(const G4Nucleus &)=default
G4ThreeVector GetFermiMomentum()
Definition G4Nucleus.cc:538
G4bool operator==(const G4Nucleus &right) const
Definition G4Nucleus.hh:67