Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QMDNucleus.cc
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// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
27//
28#include <numeric>
29
30#include "G4QMDNucleus.hh"
31#include "G4SystemOfUnits.hh"
32#include "G4Proton.hh"
33#include "G4Neutron.hh"
34#include "G4NucleiProperties.hh"
35
37{
39 hbc = parameters->Get_hbc();
40}
41
42
43
45{
46 ;
47}
48
49
51{
52 G4LorentzVector p( 0 );
53 std::vector< G4QMDParticipant* >::iterator it;
54 for ( it = participants.begin() ; it != participants.end() ; it++ )
55 p += (*it)->Get4Momentum();
56
57 return p;
58}
59
60
61
63{
64
65 G4int A = 0;
66 std::vector< G4QMDParticipant* >::iterator it;
67 for ( it = participants.begin() ; it != participants.end() ; it++ )
68 {
69 if ( (*it)->GetDefinition() == G4Proton::Proton()
70 || (*it)->GetDefinition() == G4Neutron::Neutron() )
71 A++;
72 }
73
74 return A;
75}
76
77
78
80{
81 G4int Z = 0;
82 std::vector< G4QMDParticipant* >::iterator it;
83 for ( it = participants.begin() ; it != participants.end() ; it++ )
84 {
85 if ( (*it)->GetDefinition() == G4Proton::Proton() )
86 Z++;
87 }
88 return Z;
89}
90
91
92
94{
95
97
98 if ( mass == 0.0 )
99 {
100
102 G4int A = GetMassNumber();
103 G4int N = A - Z;
104
105// Weizsacker-Bethe
106
107 G4double Av = 16*MeV;
108 G4double As = 17*MeV;
109 G4double Ac = 0.7*MeV;
110 G4double Asym = 23*MeV;
111
112 G4double BE = Av * A
113 - As * std::pow ( G4double ( A ) , 2.0/3.0 )
114 - Ac * Z*Z/std::pow ( G4double ( A ) , 1.0/3.0 )
115 - Asym * ( N - Z )* ( N - Z ) / A;
116
117 mass = Z * G4Proton::Proton()->GetPDGMass()
119 - BE;
120
121 }
122
123 return mass;
124}
125
126
127
129{
130
131 //G4cout << "CalEnergyAndAngularMomentumInCM " << this->GetAtomicNumber() << " " << GetMassNumber() << G4endl;
132
133 G4double gamma = Get4Momentum().gamma();
134 G4ThreeVector beta = Get4Momentum().v()/ Get4Momentum().e();
135
136 G4ThreeVector pcm0( 0.0 ) ;
137
139 pcm.resize( n );
140
141 for ( G4int i= 0; i < n ; i++ )
142 {
144
145 G4double trans = gamma / ( gamma + 1.0 ) * p_i * beta;
146 pcm[i] = p_i - trans*beta;
147
148 pcm0 += pcm[i];
149 }
150
151 pcm0 = pcm0 / double ( n );
152
153 //G4cout << "pcm0 " << pcm0 << G4endl;
154
155 for ( G4int i= 0; i < n ; i++ )
156 {
157 pcm[i] += -pcm0;
158 //G4cout << "pcm " << i << " " << pcm[i] << G4endl;
159 }
160
161
162 G4double tmass = 0;
163 G4ThreeVector rcm0( 0.0 ) ;
164 rcm.resize( n );
165 es.resize( n );
166
167 for ( G4int i= 0; i < n ; i++ )
168 {
170 G4double trans = gamma / ( gamma + 1.0 ) * ri * beta;
171
172 es[i] = std::sqrt ( std::pow ( GetParticipant( i )->GetMass() , 2 ) + pcm[i]*pcm[i] );
173
174 rcm[i] = ri + trans*beta;
175
176 rcm0 += rcm[i]*es[i];
177
178 tmass += es[i];
179 }
180
181 rcm0 = rcm0/tmass;
182
183 for ( G4int i= 0; i < n ; i++ )
184 {
185 rcm[i] += -rcm0;
186 //G4cout << "rcm " << i << " " << rcm[i] << G4endl;
187 }
188
189// Angluar momentum
190
191 G4ThreeVector rl ( 0.0 );
192 for ( G4int i= 0; i < n ; i++ )
193 {
194 rl += rcm[i].cross ( pcm[i] );
195 }
196
197 jj = int ( std::sqrt ( rl*rl / hbc ) + 0.5 );
198
199
200// kinetic energy per nucleon in CM
201
202 G4double totalMass = 0.0;
203 for ( G4int i= 0; i < n ; i++ )
204 {
205 // following two lines are equivalent
206 //totalMass += GetParticipant( i )->GetDefinition()->GetPDGMass()/GeV;
207 totalMass += GetParticipant( i )->GetMass();
208 }
209
210 kineticEnergyPerNucleon = ( std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass )/n;
211
212// Total (not per nucleion ) Binding Energy
213 bindingEnergy = ( std::accumulate ( es.begin() , es.end() , 0.0 ) -totalMass ) + potentialEnergy;
214
215 //G4cout << "KineticEnergyPerNucleon in GeV " << kineticEnergyPerNucleon << G4endl;
216 //G4cout << "KineticEnergySum in GeV " << std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass << G4endl;
217 //G4cout << "PotentialEnergy in GeV " << potentialEnergy << G4endl;
218 //G4cout << "BindingEnergy in GeV " << bindingEnergy << G4endl;
219 //G4cout << "G4BindingEnergy in GeV " << G4NucleiProperties::GetBindingEnergy( GetAtomicNumber() , GetMassNumber() )/GeV << G4endl;
220
221 excitationEnergy = bindingEnergy + G4NucleiProperties::GetBindingEnergy( GetMassNumber() , GetAtomicNumber() )/GeV;
222 //G4cout << "excitationEnergy in GeV " << excitationEnergy << G4endl;
223 if ( excitationEnergy < 0 ) excitationEnergy = 0.0;
224
225}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
Hep3Vector v() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetBindingEnergy(const G4int A, const G4int Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4int GetAtomicNumber()
Definition: G4QMDNucleus.cc:79
G4int GetMassNumber()
Definition: G4QMDNucleus.cc:62
void CalEnergyAndAngularMomentumInCM()
G4double GetNuclearMass()
Definition: G4QMDNucleus.cc:93
G4LorentzVector Get4Momentum()
Definition: G4QMDNucleus.cc:50
static G4QMDParameters * GetInstance()
G4double Get_hbc()
G4ThreeVector GetPosition()
G4ThreeVector GetMomentum()
G4QMDParticipant * GetParticipant(G4int i)
Definition: G4QMDSystem.hh:62
std::vector< G4QMDParticipant * > participants
Definition: G4QMDSystem.hh:72
G4int GetTotalNumberOfParticipant()
Definition: G4QMDSystem.hh:60