The Elk Code
 
Loading...
Searching...
No Matches
radnucl.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: radnucl
8! !INTERFACE:
9elemental real(8) function radnucl(z)
10! !INPUT/OUTPUT PARAMETERS:
11! z : atomic number (in,real)
12! !DESCRIPTION:
13! Computes an approximate nuclear charge radius from the atomic number $Z$.
14! The nuclear mass number, $A$, is 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)], and the nuclear charge radius can be determined from
18! $$ r=\left(r_0+\frac{r_1}{A^{2/3}}+\frac{r_2}{A^{4/3}}\right)A^{1/3}, $$
19! where $r_0=0.9071$, $r_1=1.105$ and $r_2=-0.548$ [I. Angeli, {\it Atomic
20! Data and Nuclear Data Tables} {\bf 87}, 185 (2004)].
21!
22! !REVISION HISTORY:
23! Created October 2011 (JKD)
24!EOP
25!BOC
26implicit none
27! arguments
28real(8), intent(in) :: z
29! local variables
30! coefficients for computing mass number
31real(8), parameter :: c2=4.467d-3, c1=2.163d0, c0=-1.168d0
32! coefficients for computing charge radius (fm)
33real(8), parameter :: r0=0.9071d0, r1=1.105d0, r2=-0.548d0
34! Bohr radius in SI units (CODATA 2018)
35real(8), parameter :: br_si=0.529177210903d-10
36real(8) za,a,a13,a23,a43
37za=abs(z)
38! approximate nuclear mass number
39if (za <= 1.d0) then
40 a=1.d0
41else
42 a=abs(c2*za**2+c1*za+c0)
43end if
44! approximate nuclear charge radius
45a13=a**(1.d0/3.d0)
46a23=a13**2
47a43=a13*a
48radnucl=(r0+r1/a23+r2/a43)*a13
49radnucl=radnucl*1.d-15/br_si
50end function
51!EOC
52
elemental real(8) function radnucl(z)
Definition radnucl.f90:10