Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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
34#ifndef G4Nucleus_h
35#define G4Nucleus_h 1
36// Class Description
37// This class knows how to describe a nucleus;
38// to be used in your physics implementation (not physics list) in case you need this physics.
39// Class Description - End
40
41
42#include "globals.hh"
43#include "G4ThreeVector.hh"
44#include "G4ParticleTypes.hh"
45#include "G4ReactionProduct.hh"
46#include "G4DynamicParticle.hh"
48#include "Randomize.hh"
49
51{
52 public:
53
54 G4Nucleus();
55 G4Nucleus(const G4double A, const G4double Z);
56 G4Nucleus(const G4int A, const G4int Z);
57 G4Nucleus(const G4Material* aMaterial);
58
59 ~G4Nucleus();
60
61 inline G4Nucleus( const G4Nucleus &right )
62 { *this = right; }
63
64 inline G4Nucleus& operator = (const G4Nucleus& right)
65 {
66 if (this != &right) {
67 theA=right.theA;
68 theZ=right.theZ;
69 aEff=right.aEff;
70 zEff=right.zEff;
71 fIsotope = right.fIsotope;
72 pnBlackTrackEnergy=right.pnBlackTrackEnergy;
73 dtaBlackTrackEnergy=right.dtaBlackTrackEnergy;
74 pnBlackTrackEnergyfromAnnihilation =
75 right.pnBlackTrackEnergyfromAnnihilation;
76 dtaBlackTrackEnergyfromAnnihilation =
77 right.dtaBlackTrackEnergyfromAnnihilation;
78 theTemp = right.theTemp;
79 excitationEnergy = right.excitationEnergy;
80 momentum = right.momentum;
81 fermiMomentum = right.fermiMomentum;
82 }
83 return *this;
84 }
85
86 inline G4bool operator==( const G4Nucleus &right ) const
87 { return ( this == (G4Nucleus *) &right ); }
88
89 inline G4bool operator!=( const G4Nucleus &right ) const
90 { return ( this != (G4Nucleus *) &right ); }
91
92 void ChooseParameters( const G4Material *aMaterial );
93
94 void SetParameters( const G4double A, const G4double Z );
95 void SetParameters( const G4int A, const G4int Z );
96
97/*
98#ifndef G4Hadr_Nucleus_IntegerAZ
99//deprecated Jan 2010, GF
100 inline G4double GetN() const
101 { return aEff; }
102
103 inline G4double GetZ() const
104 { return zEff; }
105#endif
106//to be replaced by new
107*/
108
109 inline G4int GetA_asInt() const
110 { return theA; }
111
112 inline G4int GetN_asInt() const
113 { return theA-theZ; }
114
115 inline G4int GetZ_asInt() const
116 { return theZ; }
117//... \GF
118
119 inline const G4Isotope* GetIsotope()
120 { return fIsotope; }
121
122 inline void SetIsotope(const G4Isotope* iso)
123 {
124 fIsotope = iso;
125 if(iso) {
126 theZ = iso->GetZ();
127 theA = iso->GetN();
128 aEff = theA;
129 zEff = theZ;
130 }
131 }
132
134
135 G4double AtomicMass( const G4double A, const G4double Z ) const;
136 G4double AtomicMass( const G4int A, const G4int Z ) const;
137
138 G4double GetThermalPz( const G4double mass, const G4double temp ) const;
139
141
143
144 G4double Cinema( G4double kineticEnergy );
145
146 G4double EvaporationEffects( G4double kineticEnergy );
147
149
151 { return pnBlackTrackEnergy; }
152
154 { return dtaBlackTrackEnergy; }
155
157 { return pnBlackTrackEnergyfromAnnihilation; }
158
160 { return dtaBlackTrackEnergyfromAnnihilation; }
161
162// ****************** methods introduced by ChV ***********************
163 // return fermi momentum
165
166/*
167 // return particle to be absorbed.
168 G4DynamicParticle* ReturnAbsorbingParticle(G4double weight);
169*/
170
171 // final nucleus fragmentation. Return List of particles
172 // which should be used for further tracking.
174
175
176 // excitation Energy...
177 void AddExcitationEnergy(G4double anEnergy);
178
179
180 // momentum of absorbed Particles ..
181 void AddMomentum(const G4ThreeVector aMomentum);
182
183 // return excitation Energy
184 G4double GetEnergyDeposit() {return excitationEnergy; }
185
186
187
188// ****************************** end ChV ******************************
189
190
191 private:
192
193 G4int theA;
194 G4int theZ;
195 G4double aEff; // effective atomic weight
196 G4double zEff; // effective atomic number
197
198 const G4Isotope* fIsotope;
199
200 G4double pnBlackTrackEnergy; // the kinetic energy available for
201 // proton/neutron black track particles
202 G4double dtaBlackTrackEnergy; // the kinetic energy available for
203 // deuteron/triton/alpha particles
204 G4double pnBlackTrackEnergyfromAnnihilation;
205 // kinetic energy available for proton/neutron black
206 // track particles based on baryon annihilation
207 G4double dtaBlackTrackEnergyfromAnnihilation;
208 // kinetic energy available for deuteron/triton/alpha
209 // black track particles based on baryon annihilation
210
211
212// ************************** member variables by ChV *******************
213 // Excitation Energy leading to evaporation or deexcitation.
214 G4double excitationEnergy;
215
216 // Momentum, accumulated by absorbing Particles
217 G4ThreeVector momentum;
218
219 // Fermi Gas model: at present, we assume constant nucleon density for all
220 // nuclei. The radius of a nucleon is taken to be 1 fm.
221 // see for example S.Fl"ugge, Encyclopedia of Physics, Vol XXXIX,
222 // Structure of Atomic Nuclei (Berlin-Gottingen-Heidelberg, 1957) page 426.
223
224 // maximum momentum possible from fermi gas model:
225 G4double fermiMomentum;
226 G4double theTemp; // temperature
227// ****************************** end ChV ******************************
228
229 };
230
231#endif
232
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4int GetZ() const
Definition: G4Isotope.hh:91
G4int GetN() const
Definition: G4Isotope.hh:94
void AddExcitationEnergy(G4double anEnergy)
Definition: G4Nucleus.cc:435
G4Nucleus(const G4Nucleus &right)
Definition: G4Nucleus.hh:61
G4double GetThermalPz(const G4double mass, const G4double temp) const
Definition: G4Nucleus.cc:254
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
void ChooseParameters(const G4Material *aMaterial)
Definition: G4Nucleus.cc:158
G4double GetAnnihilationPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:156
G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg)
Definition: G4Nucleus.cc:323
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:119
G4double GetPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:150
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4double GetAnnihilationDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:159
G4bool operator!=(const G4Nucleus &right) const
Definition: G4Nucleus.hh:89
G4int GetN_asInt() const
Definition: G4Nucleus.hh:112
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:227
G4double GetEnergyDeposit()
Definition: G4Nucleus.hh:184
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
G4ReactionProductVector * Fragmentate()
Definition: G4Nucleus.cc:424
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
G4double GetDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:153
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:240
void AddMomentum(const G4ThreeVector aMomentum)
Definition: G4Nucleus.cc:430
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
G4Nucleus & operator=(const G4Nucleus &right)
Definition: G4Nucleus.hh:64
G4ThreeVector GetFermiMomentum()
Definition: G4Nucleus.cc:398
G4bool operator==(const G4Nucleus &right) const
Definition: G4Nucleus.hh:86