The Elk Code
 
Loading...
Searching...
No Matches
rdiracint.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2015 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: rdiracint
8! !INTERFACE:
9pure subroutine rdiracint(sol,kpa,e,nr,r,vr,nn,g0,g1,f0,f1)
10! !INPUT/OUTPUT PARAMETERS:
11! sol : speed of light in atomic units (in,real)
12! kpa : quantum number kappa (in,integer)
13! e : energy (in,real)
14! nr : number of radial mesh points (in,integer)
15! r : radial mesh (in,real(nr))
16! vr : potential on radial mesh (in,real(nr))
17! nn : number of nodes (out,integer)
18! g0 : major component multiplied by r (out,real(nr))
19! g1 : radial derivative of g0 (out,real(nr))
20! f0 : minor component multiplied by r (out,real(nr))
21! f1 : radial derivative of f0 (out,real(nr))
22! !DESCRIPTION:
23! Integrates the radial Dirac equation from $r=0$ outwards. This involves
24! using the predictor-corrector method to solve the coupled first-order
25! equations (in atomic units)
26! \begin{align*}
27! \left(\frac{d}{dr}+\frac{\kappa}{r}\right)G_{\kappa}&=\frac{1}{c}
28! \{2E_0+E-V\}F_{\kappa}\\
29! \left(\frac{d}{dr}-\frac{\kappa}{r}\right)F_{\kappa}&=
30! -\frac{1}{c}\{E-V\}G_{\kappa},
31! \end{align*}
32! where $G_{\kappa}=rg_{\kappa}$ and $F_{\kappa}=rf_{\kappa}$ are the major
33! and minor components multiplied by $r$, respectively; $V$ is the external
34! potential; $E_0$ is the electron rest energy; $E$ is the eigen energy
35! (excluding $E_0$); and $\kappa=l$ for $j=l-\frac{1}{2}$ or $\kappa=-(l+1)$
36! for $j=l+\frac{1}{2}$.
37!
38! !REVISION HISTORY:
39! Created September 2002 (JKD)
40! Polynomial order fixed to 3, September 2013 (JKD)
41!EOP
42!BOC
43implicit none
44! arguments
45real(8), intent(in) :: sol
46integer, intent(in) :: kpa
47real(8), intent(in) :: e
48integer, intent(in) :: nr
49real(8), intent(in) :: r(nr),vr(nr)
50integer, intent(out) :: nn
51real(8), intent(out) :: g0(nr),g1(nr),f0(nr),f1(nr)
52! local variables
53integer ir,ir0,i
54real(8) ci,e0,t1,t2,t3,t4
55! inverse speed of light
56ci=1.d0/sol
57! electron rest energy
58e0=sol**2
59t1=2.d0*e0+e
60! determine the r -> 0 boundary values of F and G
61t2=dble(kpa)/r(1)
62t3=ci*(t1-vr(1))
63t4=ci*(vr(1)-e)
64f0(1)=1.d0
65g0(1)=-t2/t4
66g1(1)=t3-t2*g0(1)
67! extrapolate to the first four points
68g1(2:4)=g1(1)
69f1(1:4)=0.d0
70nn=0
71do ir=2,nr
72 t2=dble(kpa)/r(ir)
73 t3=ci*(t1-vr(ir))
74 t4=ci*(vr(ir)-e)
75 ir0=max(ir-3,1)
76! estimate the derivatives
77 g1(ir)=poly3(r(ir0),g1(ir0),r(ir))
78 f1(ir)=poly3(r(ir0),f1(ir0),r(ir))
79 do i=1,8
80! integrate to find wavefunction
81 g0(ir)=poly4i(r(ir0),g1(ir0),r(ir))+g0(ir0)
82 f0(ir)=poly4i(r(ir0),f1(ir0),r(ir))+f0(ir0)
83! compute the derivatives
84 g1(ir)=t3*f0(ir)-t2*g0(ir)
85 f1(ir)=t4*g0(ir)+t2*f0(ir)
86 end do
87! check for overflow
88 if (abs(g0(ir)) > 1.d100) then
89! set the remaining points and return
90 g0(ir+1:nr)=g0(ir); g1(ir+1:nr)=g1(ir)
91 f0(ir+1:nr)=f0(ir); f1(ir+1:nr)=f1(ir)
92 return
93 end if
94! check for node
95 if (g0(ir-1)*g0(ir) < 0.d0) nn=nn+1
96end do
97return
98
99contains
100
101pure real(8) function poly3(xa,ya,x)
102implicit none
103! arguments
104real(8), intent(in) :: xa(3),ya(3),x
105! local variables
106real(8) x0,x1,x2,y0,y1,y2
107real(8) c1,c2,t0,t1,t2
108! evaluate the polynomial coefficients
109x0=xa(1)
110x1=xa(2)-x0; x2=xa(3)-x0
111y0=ya(1)
112y1=ya(2)-y0; y2=ya(3)-y0
113t0=1.d0/(x1*x2*(x2-x1))
114t1=x1*y2; t2=x2*y1
115c1=x2*t2-x1*t1
116c2=t1-t2
117t1=x-x0
118! evaluate the polynomial
119poly3=y0+t0*t1*(c1+c2*t1)
120end function
121
122pure real(8) function poly4i(xa,ya,x)
123implicit none
124! arguments
125real(8), intent(in) :: xa(4),ya(4),x
126! local variables
127real(8) x0,x1,x2,x3,y0,y1,y2,y3
128real(8) c1,c2,c3,t0,t1,t2,t3,t4,t5,t6
129! evaluate the polynomial coefficients
130x0=xa(1)
131x1=xa(2)-x0; x2=xa(3)-x0; x3=xa(4)-x0
132y0=ya(1)
133y1=ya(2)-y0; y2=ya(3)-y0; y3=ya(4)-y0
134t4=x1-x2; t5=x1-x3; t6=x2-x3
135t1=x1*x2*y3; t2=x2*x3*y1; t3=x1*x3
136t0=1.d0/(x2*t3*t4*t5*t6)
137t3=t3*y2
138c3=t1*t4+t2*t6-t3*t5
139t4=x1**2; t5=x2**2; t6=x3**2
140c2=t1*(t5-t4)+t2*(t6-t5)+t3*(t4-t6)
141c1=t1*(x2*t4-x1*t5)+t2*(x3*t5-x2*t6)+t3*(x1*t6-x3*t4)
142t1=x-x0
143! integrate the polynomial
144poly4i=t1*(y0+t0*t1*(0.5d0*c1+t1*(0.3333333333333333333d0*c2+0.25d0*c3*t1)))
145end function
146
147end subroutine
148!EOC
149
pure subroutine rdiracint(sol, kpa, e, nr, r, vr, nn, g0, g1, f0, f1)
Definition rdiracint.f90:10
pure real(8) function poly4i(xa, ya, x)
pure real(8) function poly3(xa, ya, x)