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