Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LightIonQMDNucleus Class Reference

#include <G4LightIonQMDNucleus.hh>

Inheritance diagram for G4LightIonQMDNucleus:

Public Member Functions

 G4LightIonQMDNucleus ()
G4LorentzVector Get4Momentum ()
G4int GetMassNumber ()
G4int GetAtomicNumber ()
void CalEnergyAndAngularMomentumInCM ()
G4double GetNuclearMass ()
void SetTotalPotential (G4double x)
G4double GetExcitationEnergy ()
G4int GetAngularMomentum ()
Public Member Functions inherited from G4QMDSystem
 G4QMDSystem ()
virtual ~G4QMDSystem ()
void SetParticipant (G4QMDParticipant *particle)
void SetSystem (G4QMDSystem *, G4ThreeVector, G4ThreeVector)
void SubtractSystem (G4QMDSystem *)
G4QMDParticipantEraseParticipant (G4int i)
void DeleteParticipant (G4int i)
void InsertParticipant (G4QMDParticipant *particle, G4int j)
G4int GetTotalNumberOfParticipant ()
G4QMDParticipantGetParticipant (G4int i)
void IncrementCollisionCounter ()
G4int GetNOCollision ()
void ShowParticipants ()
void Clear ()

Additional Inherited Members

Protected Attributes inherited from G4QMDSystem
std::vector< G4QMDParticipant * > participants

Detailed Description

Definition at line 44 of file G4LightIonQMDNucleus.hh.

Constructor & Destructor Documentation

◆ G4LightIonQMDNucleus()

G4LightIonQMDNucleus::G4LightIonQMDNucleus ( )

Definition at line 46 of file G4LightIonQMDNucleus.cc.

47{
48 G4LightIonQMDParameters* parameters = G4LightIonQMDParameters::GetInstance();
49 hbc = parameters->Get_hbc();
50
51 jj = 0; // will be calcualted in CalEnergyAndAngularMomentumInCM;
52 potentialEnergy = 0.0; // will be set through set method
53 excitationEnergy = 0.0;
54
55 // Following Parameters are added (20230309)
56 wl = parameters->Get_wl();
57 cl = parameters->Get_cl();
58 rho0 = parameters->Get_rho0();
59 gamm = parameters->Get_gamm();
60 eta = parameters->Get_eta(); // Skyrme-QMD
61 kappas = parameters->Get_kappas(); // Skyrme-QMD
62
63 cpw = parameters->Get_cpw();
64 cph = parameters->Get_cph();
65 cpc = parameters->Get_cpc();
66
67 c0 = parameters->Get_c0();
68 c3 = parameters->Get_c3();
69 cs = parameters->Get_cs();
70 g0 = parameters->Get_g0(); // Skyrme-QMD
71 g0iso = parameters->Get_g0iso(); // Skyrme-QMD
72 gtau0 = parameters->Get_gtau0(); // Skyrme-QMD
73
74 // distance
75 c0w = 1.0/4.0/wl;
76 //c3w = 1.0/4.0/wl; //no need
77 c0sw = std::sqrt( c0w );
78 clw = 2.0 / std::sqrt ( 4.0 * pi * wl );
79
80 // graduate
81 c0g = - c0 / ( 2.0 * wl );
82 c3g = - c3 / ( 4.0 * wl ) * gamm;
83 csg = - cs / ( 2.0 * wl );
84 pag = gamm - 1;
85 pag_tau = eta - 1; // Skyrme-QMD
86 cg0 = - g0 / ( 2.0 * wl ); // Skyrme-QMD
87 cgtau0 = - gtau0 / ( 4.0 * wl ) * eta; // Skyrme-QMD
88}
static G4LightIonQMDParameters * GetInstance()

Member Function Documentation

◆ CalEnergyAndAngularMomentumInCM()

void G4LightIonQMDNucleus::CalEnergyAndAngularMomentumInCM ( )

Definition at line 180 of file G4LightIonQMDNucleus.cc.

181{
182
183 //G4cout << "CalEnergyAndAngularMomentumInCM " << this->GetAtomicNumber() << " " << GetMassNumber() << G4endl;
184
185 G4double gamma = Get4Momentum().gamma();
186 G4ThreeVector beta = Get4Momentum().v()/ Get4Momentum().e();
187
188 G4ThreeVector pcm0( 0.0 ) ;
189
191 pcm.resize( n );
192
193 for ( G4int i= 0; i < n ; i++ )
194 {
196
197 G4double trans = gamma / ( gamma + 1.0 ) * p_i * beta;
198 pcm[i] = p_i - trans*beta;
199
200 pcm0 += pcm[i];
201 }
202
203 pcm0 = pcm0 / double ( n );
204
205 //G4cout << "pcm0 " << pcm0 << G4endl;
206
207 for ( G4int i= 0; i < n ; i++ )
208 {
209 pcm[i] += -pcm0;
210 //G4cout << "pcm " << i << " " << pcm[i] << G4endl;
211 }
212
213
214 G4double tmass = 0;
215 G4ThreeVector rcm0( 0.0 ) ;
216 rcm.resize( n );
217 es.resize( n );
218
219 // binding energy should be evaluated with a relativistic version: 20230308 by Y-H. Sato and A. Haga
220 for ( G4int i= 0; i < n ; i++ )
221 {
223 G4double trans = gamma / ( gamma + 1.0 ) * ri * beta;
224 G4double nucpote = GetNuclPotential( i );
225
226 es[i] = std::sqrt ( G4Pow::GetInstance()->powN ( GetParticipant( i )->GetMass() , 2 ) + pcm[i]*pcm[i] + 2.0*GetParticipant( i )->GetMass()*nucpote) - GetParticipant( i )->GetMass(); //R-JQMD
227
228 rcm[i] = ri + trans*beta;
229
230 rcm0 += rcm[i]*es[i];
231
232 tmass += es[i];
233 }
234
235 rcm0 = rcm0/tmass;
236
237 for ( G4int i= 0; i < n ; i++ )
238 {
239 rcm[i] += -rcm0;
240 //G4cout << "rcm " << i << " " << rcm[i] << G4endl;
241 }
242
243// Angular momentum
244
245 G4ThreeVector rl ( 0.0 );
246 for ( G4int i= 0; i < n ; i++ )
247 {
248 rl += rcm[i].cross ( pcm[i] );
249 }
250
251// DHW: move hbc outside of sqrt to get correct units
252// jj = int ( std::sqrt ( rl*rl / hbc ) + 0.5 );
253
254 jj = int (std::sqrt(rl*rl)/hbc + 0.5);
255
256// kinetic energy per nucleon in CM
257
258 /*
259 G4double totalMass = 0.0;
260 for ( G4int i= 0; i < n ; i++ )
261 {
262 // following two lines are equivalent
263 //totalMass += GetParticipant( i )->GetDefinition()->GetPDGMass()/GeV;
264 totalMass += GetParticipant( i )->GetMass();
265 }
266 */
267
268 //G4double kineticEnergyPerNucleon = ( std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass )/n;
269
270// Total (not per nucleion ) Binding Energy
271 // relativistic version Y-H. Sato and A. Haga 20230309
272 G4double bindingEnergy = ( std::accumulate ( es.begin() , es.end() , 0.0 ) );
273
274 //G4cout << "n " << n << "totalpote " << totalpote << " " << potentialEnergy << " " << bindingEnergy << G4endl;
275 //G4cout << "KineticEnergyPerNucleon in GeV " << kineticEnergyPerNucleon << G4endl;
276 //G4cout << "KineticEnergySum in GeV " << std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass << G4endl;
277 //G4cout << "PotentialEnergy in GeV " << potentialEnergy << G4endl;
278 //G4cout << "BindingEnergy in GeV " << bindingEnergy << G4endl;
279 //G4cout << "G4BindingEnergy in GeV " << G4NucleiProperties::GetBindingEnergy( GetAtomicNumber() , GetMassNumber() )/GeV << G4endl;
280
282 if ( excitationEnergy < 0 ) excitationEnergy = 0.0;
283
284 }
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
Hep3Vector v() const
G4LorentzVector Get4Momentum()
static G4double GetBindingEnergy(const G4int A, const G4int Z)
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4ThreeVector GetPosition()
G4ThreeVector GetMomentum()
G4QMDParticipant * GetParticipant(G4int i)
G4int GetTotalNumberOfParticipant()
G4double bindingEnergy(G4int A, G4int Z)

Referenced by G4LightIonQMDReaction::ApplyYourself(), and G4LightIonQMDMeanField::SetNucleus().

◆ Get4Momentum()

G4LorentzVector G4LightIonQMDNucleus::Get4Momentum ( )

Definition at line 98 of file G4LightIonQMDNucleus.cc.

99{
100 G4LorentzVector p( 0 );
101 std::vector< G4QMDParticipant* >::iterator it;
102 for ( it = participants.begin() ; it != participants.end() ; it++ )
103 p += (*it)->Get4Momentum();
104
105 return p;
106}
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4QMDParticipant * > participants

Referenced by CalEnergyAndAngularMomentumInCM().

◆ GetAngularMomentum()

G4int G4LightIonQMDNucleus::GetAngularMomentum ( )
inline

Definition at line 66 of file G4LightIonQMDNucleus.hh.

66{ return jj; };

◆ GetAtomicNumber()

G4int G4LightIonQMDNucleus::GetAtomicNumber ( )

Definition at line 131 of file G4LightIonQMDNucleus.cc.

132{
133 G4int Z = 0;
134 std::vector< G4QMDParticipant* >::iterator it;
135 for ( it = participants.begin() ; it != participants.end() ; it++ )
136 {
137 if ( (*it)->GetDefinition() == G4Proton::Proton() )
138 Z++;
139 }
140 return Z;
141}
static G4Proton * Proton()
Definition G4Proton.cc:90

Referenced by CalEnergyAndAngularMomentumInCM(), and GetNuclearMass().

◆ GetExcitationEnergy()

G4double G4LightIonQMDNucleus::GetExcitationEnergy ( )
inline

Definition at line 64 of file G4LightIonQMDNucleus.hh.

64{ return excitationEnergy; };

◆ GetMassNumber()

G4int G4LightIonQMDNucleus::GetMassNumber ( )

Definition at line 110 of file G4LightIonQMDNucleus.cc.

111{
112
113 G4int A = 0;
114 std::vector< G4QMDParticipant* >::iterator it;
115 for ( it = participants.begin() ; it != participants.end() ; it++ )
116 {
117 if ( (*it)->GetDefinition() == G4Proton::Proton()
118 || (*it)->GetDefinition() == G4Neutron::Neutron() )
119 A++;
120 }
121
122 if ( A == 0 ) {
123 throw G4HadronicException(__FILE__, __LINE__, "G4LightIonQMDNucleus has the mass number of 0!");
124 }
125
126 return A;
127}
const G4double A[17]
static G4Neutron * Neutron()
Definition G4Neutron.cc:101

Referenced by CalEnergyAndAngularMomentumInCM(), G4LightIonQMDGroundStateNucleus::G4LightIonQMDGroundStateNucleus(), and GetNuclearMass().

◆ GetNuclearMass()

G4double G4LightIonQMDNucleus::GetNuclearMass ( )

Definition at line 145 of file G4LightIonQMDNucleus.cc.

146{
147
149
150 if ( mass == 0.0 )
151 {
152
155 G4int N = A - Z;
156
157// Weizsacker-Bethe
158
159 G4double Av = 16*MeV;
160 G4double As = 17*MeV;
161 G4double Ac = 0.7*MeV;
162 G4double Asym = 23*MeV;
163
164 G4double BE = Av * A
165 - As * G4Pow::GetInstance()->A23 ( G4double ( A ) )
166 - Ac * Z*Z/G4Pow::GetInstance()->A13 ( G4double ( A ) )
167 - Asym * ( N - Z )* ( N - Z ) / A;
168
169 mass = Z * G4Proton::Proton()->GetPDGMass()
171 - BE;
172
173 }
174
175 return mass;
176}
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double A13(G4double A) const
Definition G4Pow.cc:116
G4double A23(G4double A) const
Definition G4Pow.hh:131
#define N
Definition crc32.c:57

◆ SetTotalPotential()

void G4LightIonQMDNucleus::SetTotalPotential ( G4double x)
inline

Definition at line 63 of file G4LightIonQMDNucleus.hh.

63{ potentialEnergy = x; };

Referenced by G4LightIonQMDReaction::ApplyYourself(), and G4LightIonQMDMeanField::SetNucleus().


The documentation for this class was generated from the following files: