source: trunk/source/processes/hadronic/models/qmd/src/G4QMDNucleus.cc @ 1055

Last change on this file since 1055 was 962, checked in by garnier, 15 years ago

update processes

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