The Elk Code
massnucl.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: massnucl
8 ! !INTERFACE:
9 elemental real(8) function massnucl(z)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! z : atomic number (in,real)
12 ! !DESCRIPTION:
13 ! Computes an approximate nuclear mass from the atomic number $Z$. The nuclear
14 ! mass number, $A$, is first estimated using
15 ! $$ A=4.467\times 10^{-3}Z^2+2.163 Z-1.168, $$
16 ! [D. Andrae in {\it Relativistic Electronic Structure Theory - Fundamentals}
17 ! {\bf 11}, 203 (2002)]. Then the nuclear mass can be determined from:
18 ! $$ M=Z m_p+N m_n-\frac{B}{c^2}, $$
19 ! where $m_p$ is the proton mass, $m_n$ is the neutron mass and $B$ is the
20 ! nuclear binding energy. The binding energy is approximated by the
21 ! Weizs\"{a}cker formula:
22 ! $$ B=a_V A-a_S A^{2/3}-a_C Z^2 A^{-1/3}-a_{\rm sym}(Z-N)^2A^{-1}
23 ! +B_p+B_{\rm shell}. $$
24 ! See F. Yang and J. H. Hamilton in {\it Modern Atomic and Nuclear Physics},
25 ! Revised Edition 2010, for details on the quantities in this formula. In this
26 ! implementation, $B_p$ and $B_{\rm shell}$ are set to zero.
27 !
28 ! !REVISION HISTORY:
29 ! Created February 2014 (JKD)
30 !EOP
31 !BOC
32 implicit none
33 ! arguments
34 real(8), intent(in) :: z
35 ! local variables
36 ! coefficients for computing mass number
37 real(8), parameter :: c2=4.467d-3, c1=2.163d0, c0=-1.168d0
38 ! Weizsacker coefficients in MeV
39 real(8), parameter :: av=15.8d0, as=18.3d0, ac=0.72d0, asym=23.2d0
40 ! proton and neutron masses in MeV/c² (CODATA 2018)
41 real(8), parameter :: mp=938.27208816d0
42 real(8), parameter :: mn=939.56542052d0
43 ! atomic mass unit in MeV/c² (CODATA 2018)
44 real(8), parameter :: mu=931.49410242d0
45 real(8) za,n,a,b
46 za=abs(z)
47 ! approximate nuclear mass number
48 if (za <= 1.d0) then
49  a=1.d0
50 else
51  a=abs(c2*za**2+c1*za+c0)
52 end if
53 n=a-za
54 b=av*a-as*a**(2.d0/3.d0)-ac*(za**2)/a**(1.d0/3.d0)-asym*(za-n)**2/a
55 massnucl=(za*mp+n*mn-b)/mu
56 end function
57 !EOC
58 
elemental real(8) function massnucl(z)
Definition: massnucl.f90:10