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