The Elk Code
Loading...
Searching...
No Matches
potnucl.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: potnucl
8
! !INTERFACE:
9
subroutine
potnucl
(ptnucl,nr,r,zn,vn)
10
! !INPUT/OUTPUT PARAMETERS:
11
! ptnucl : .true. if the nucleus is a point charge (in,logical)
12
! nr : number of radial mesh points (in,integer)
13
! r : radial mesh (in,real(nr))
14
! zn : nuclear charge (in,real)
15
! vn : potential on radial mesh (out,real(nr))
16
! !DESCRIPTION:
17
! Computes the nuclear Coulomb potential on a radial mesh. The nuclear radius
18
! $R$ is estimated from the nuclear charge $Z$ and the potential is given by
19
! $$ V(r)=\begin{cases}
20
! Z(3R^2-r^2)/2R^3 & r<R \\
21
! Z/r & r\ge R\end{cases} $$
22
! assuming that the nucleus is a uniformly charged sphere. If {\tt ptnucl} is
23
! {\tt .true.} then the nucleus is treated as a point particle.
24
!
25
! !REVISION HISTORY:
26
! Created January 2009 (JKD)
27
!EOP
28
!BOC
29
implicit none
30
! arguments
31
logical
,
intent(in)
:: ptnucl
32
integer
,
intent(in)
:: nr
33
real
(8),
intent(in)
:: r(nr),zn
34
real
(8),
intent(out)
:: vn(nr)
35
! local variables
36
real
(8) rn,t1,t2
37
! external functions
38
real
(8),
external
:: radnucl
39
if
(zn == 0.d0)
then
40
vn(1:nr)=0.d0
41
return
42
end if
43
if
(ptnucl)
then
44
! nucleus is taken to be a point particle
45
vn(1:nr)=zn/r(1:nr)
46
else
47
! approximate nuclear radius
48
rn=
radnucl
(zn)
49
t1=zn/(2.d0*rn**3)
50
t2=3.d0*rn**2
51
where
(r(:) < rn)
52
vn(:)=t1*(t2-r(:)**2)
53
elsewhere
54
vn(:)=zn/r(:)
55
end where
56
end if
57
end subroutine
58
!EOC
59
potnucl
subroutine potnucl(ptnucl, nr, r, zn, vn)
Definition
potnucl.f90:10
radnucl
elemental real(8) function radnucl(z)
Definition
radnucl.f90:10
potnucl.f90
Generated by
1.9.8