Geant4-11
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 inline G4Nucleus( const G4Nucleus &right )
63 { *this = right; }
64
65 inline G4Nucleus& operator = (const G4Nucleus& right)
66 {
67 if (this != &right) {
68 theA=right.theA;
69 theZ=right.theZ;
70 theL=right.theL;
71 aEff=right.aEff;
72 zEff=right.zEff;
73 fIsotope = right.fIsotope;
80 theTemp = right.theTemp;
82 momentum = right.momentum;
84 }
85 return *this;
86 }
87
88 inline G4bool operator==( const G4Nucleus &right ) const
89 { return ( this == (G4Nucleus *) &right ); }
90
91 inline G4bool operator!=( const G4Nucleus &right ) const
92 { return ( this != (G4Nucleus *) &right ); }
93
94 void ChooseParameters( const G4Material *aMaterial );
95
96 void SetParameters( const G4double A, const G4double Z, const G4int numberOfLambdas = 0 );
97 void SetParameters( const G4int A, const G4int Z, const G4int numberOfLambdas = 0 );
98
99 inline G4int GetA_asInt() const
100 { return theA; }
101
102 inline G4int GetN_asInt() const
103 { return theA-theZ; }
104
105 inline G4int GetZ_asInt() const
106 { return theZ; }
107
108 inline G4int GetL() const // Number of Lambdas (in the case of a hypernucleus)
109 { return theL; }
110
111 inline const G4Isotope* GetIsotope()
112 { return fIsotope; }
113
114 inline void SetIsotope(const G4Isotope* iso)
115 {
116 fIsotope = iso;
117 if(iso) {
118 theZ = iso->GetZ();
119 theA = iso->GetN();
120 theL = 0;
121 aEff = theA;
122 zEff = theZ;
123 }
124 }
125
127
128 G4double AtomicMass( const G4double A, const G4double Z, const G4int numberOfLambdas = 0 ) const;
129 G4double AtomicMass( const G4int A, const G4int Z, const G4int numberOfLambdas = 0 ) const;
130
131 G4double GetThermalPz( const G4double mass, const G4double temp ) const;
132
134
136
137 G4double Cinema( G4double kineticEnergy );
138
139 G4double EvaporationEffects( G4double kineticEnergy );
140
142
144 { return pnBlackTrackEnergy; }
145
147 { return dtaBlackTrackEnergy; }
148
151
154
155// ****************** methods introduced by ChV ***********************
156 // return fermi momentum
158
159/*
160 // return particle to be absorbed.
161 G4DynamicParticle* ReturnAbsorbingParticle(G4double weight);
162*/
163
164 // final nucleus fragmentation. Return List of particles
165 // which should be used for further tracking.
167
168
169 // excitation Energy...
170 void AddExcitationEnergy(G4double anEnergy);
171
172
173 // momentum of absorbed Particles ..
174 void AddMomentum(const G4ThreeVector aMomentum);
175
176 // return excitation Energy
178
179
180
181// ****************************** end ChV ******************************
182
183
184 private:
185
188 G4int theL; // Number of Lambdas (in the case of hypernucleus)
189 G4double aEff; // effective atomic weight
190 G4double zEff; // effective atomic number
191
193
194 G4double pnBlackTrackEnergy; // the kinetic energy available for
195 // proton/neutron black track particles
196 G4double dtaBlackTrackEnergy; // the kinetic energy available for
197 // deuteron/triton/alpha particles
199 // kinetic energy available for proton/neutron black
200 // track particles based on baryon annihilation
202 // kinetic energy available for deuteron/triton/alpha
203 // black track particles based on baryon annihilation
204
205
206// ************************** member variables by ChV *******************
207 // Excitation Energy leading to evaporation or deexcitation.
209
210 // Momentum, accumulated by absorbing Particles
212
213 // Fermi Gas model: at present, we assume constant nucleon density for all
214 // nuclei. The radius of a nucleon is taken to be 1 fm.
215 // see for example S.Fl"ugge, Encyclopedia of Physics, Vol XXXIX,
216 // Structure of Atomic Nuclei (Berlin-Gottingen-Heidelberg, 1957) page 426.
217
218 // maximum momentum possible from fermi gas model:
220 G4double theTemp; // temperature
221// ****************************** end ChV ******************************
222
223 };
224
225#endif
226
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4int GetZ() const
Definition: G4Isotope.hh:90
G4int GetN() const
Definition: G4Isotope.hh:93
void AddExcitationEnergy(G4double anEnergy)
Definition: G4Nucleus.cc:566
G4Nucleus(const G4Nucleus &right)
Definition: G4Nucleus.hh:62
G4double GetThermalPz(const G4double mass, const G4double temp) const
Definition: G4Nucleus.cc:381
G4double aEff
Definition: G4Nucleus.hh:189
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4double dtaBlackTrackEnergy
Definition: G4Nucleus.hh:196
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:392
G4double zEff
Definition: G4Nucleus.hh:190
void ChooseParameters(const G4Material *aMaterial)
Definition: G4Nucleus.cc:265
G4double GetAnnihilationPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:149
G4int GetL() const
Definition: G4Nucleus.hh:108
G4double pnBlackTrackEnergyfromAnnihilation
Definition: G4Nucleus.hh:198
G4double AtomicMass(const G4double A, const G4double Z, const G4int numberOfLambdas=0) const
Definition: G4Nucleus.cc:357
G4double excitationEnergy
Definition: G4Nucleus.hh:208
G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg)
Definition: G4Nucleus.cc:453
const G4Isotope * GetIsotope()
Definition: G4Nucleus.hh:111
G4double GetPNBlackTrackEnergy() const
Definition: G4Nucleus.hh:143
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition: G4Nucleus.cc:307
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:499
G4double GetAnnihilationDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:152
G4bool operator!=(const G4Nucleus &right) const
Definition: G4Nucleus.hh:91
G4int GetN_asInt() const
Definition: G4Nucleus.hh:102
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:340
G4ThreeVector momentum
Definition: G4Nucleus.hh:211
G4double theTemp
Definition: G4Nucleus.hh:220
G4int theZ
Definition: G4Nucleus.hh:187
G4int theA
Definition: G4Nucleus.hh:186
G4double GetEnergyDeposit()
Definition: G4Nucleus.hh:177
G4double fermiMomentum
Definition: G4Nucleus.hh:219
G4ReactionProductVector * Fragmentate()
Definition: G4Nucleus.cc:553
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:114
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:118
G4double GetDTABlackTrackEnergy() const
Definition: G4Nucleus.hh:146
G4double dtaBlackTrackEnergyfromAnnihilation
Definition: G4Nucleus.hh:201
G4double pnBlackTrackEnergy
Definition: G4Nucleus.hh:194
void AddMomentum(const G4ThreeVector aMomentum)
Definition: G4Nucleus.cc:560
const G4Isotope * fIsotope
Definition: G4Nucleus.hh:192
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:236
G4Nucleus & operator=(const G4Nucleus &right)
Definition: G4Nucleus.hh:65
G4ThreeVector GetFermiMomentum()
Definition: G4Nucleus.cc:526
G4int theL
Definition: G4Nucleus.hh:188
G4bool operator==(const G4Nucleus &right) const
Definition: G4Nucleus.hh:88