The Elk Code
rschrodint.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: rschrodint
8 ! !INTERFACE:
9 pure subroutine rschrodint(sol,l,e,nr,r,vr,nn,p0,p1,q0,q1)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! sol : speed of light in atomic units (in,real)
12 ! l : angular momentum quantum number (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 ! p0 : radial function P (out,real(nr))
19 ! p1 : radial derivative of P (out,real(nr))
20 ! q0 : radial function Q (out,real(nr))
21 ! q1 : radial derivative of Q (out,real(nr))
22 ! !DESCRIPTION:
23 ! Integrates the scalar relativistic radial Schr\"{o}dinger equation from
24 ! $r=0$ outwards. This involves using the predictor-corrector method to solve
25 ! the coupled first-order equations (in atomic units)
26 ! \begin{align*}
27 ! \frac{d}{dr}P_l&=2MQ_l+\frac{1}{r}P_l\\
28 ! \frac{d}{dr}Q_l&=-\frac{1}{r}Q_l+\left[\frac{l(l+1)}{2Mr^2}
29 ! +(V-E)\right]P_l,
30 ! \end{align*}
31 ! where $V$ is the external potential, $E$ is the eigen energy and
32 ! $M=1+(E-V)/2c^2$. Following the convention of Koelling and Harmon,
33 ! {\it J. Phys. C: Solid State Phys.} {\bf 10}, 3107 (1977), the functions
34 ! $P_l$ and $Q_l$ are defined by
35 ! \begin{align*}
36 ! P_l&=rg_l\\
37 ! Q_l&=\frac{r}{2M}\frac{dg_l}{dr},
38 ! \end{align*}
39 ! where $g_l$ is the major component of the Dirac equation (see the routine
40 ! {\tt rdiracint}).
41 !
42 ! !REVISION HISTORY:
43 ! Created October 2003 (JKD)
44 !EOP
45 !BOC
46 implicit none
47 ! arguments
48 real(8), intent(in) :: sol
49 integer, intent(in) :: l
50 real(8), intent(in) :: e
51 integer, intent(in) :: nr
52 real(8), intent(in) :: r(nr),vr(nr)
53 integer, intent(out) :: nn
54 real(8), intent(out) :: p0(nr),p1(nr),q0(nr),q1(nr)
55 ! local variables
56 integer ir,ir0,i
57 real(8) ri,t1,t2,t3,t4
58 t1=1.d0/sol**2
59 t2=dble(l*(l+1))
60 ! determine the r -> 0 boundary values of P and Q
61 ri=1.d0/r(1)
62 t3=2.d0+t1*(e-vr(1))
63 t4=(t2*ri**2)/t3+vr(1)-e
64 q0(1)=1.d0
65 p0(1)=ri/t4
66 p1(1)=t3+p0(1)*ri
67 ! extrapolate to the first four points
68 p1(2:4)=p1(1)
69 q1(1:4)=0.d0
70 nn=0
71 do ir=2,nr
72  ri=1.d0/r(ir)
73  t3=2.d0+t1*(e-vr(ir))
74  t4=(t2*ri**2)/t3+vr(ir)-e
75  ir0=max(ir-3,1)
76 ! estimate the derivatives
77  p1(ir)=poly3(r(ir0),p1(ir0),r(ir))
78  q1(ir)=poly3(r(ir0),q1(ir0),r(ir))
79  do i=1,8
80 ! integrate to find wavefunction
81  p0(ir)=poly4i(r(ir0),p1(ir0),r(ir))+p0(ir0)
82  q0(ir)=poly4i(r(ir0),q1(ir0),r(ir))+q0(ir0)
83 ! compute the derivatives
84  p1(ir)=t3*q0(ir)+p0(ir)*ri
85  q1(ir)=t4*p0(ir)-q0(ir)*ri
86  end do
87 ! check for overflow
88  if (abs(p0(ir)) > 1.d100) then
89  p0(ir+1:nr)=p0(ir); p1(ir+1:nr)=p1(ir)
90  q0(ir+1:nr)=q0(ir); q1(ir+1:nr)=q1(ir)
91  return
92  end if
93 ! check for node
94  if (p0(ir-1)*p0(ir) < 0.d0) nn=nn+1
95 end do
96 
97 contains
98 
99 pure real(8) function poly3(xa,ya,x)
100 implicit none
101 ! arguments
102 real(8), intent(in) :: xa(3),ya(3),x
103 ! local variables
104 real(8) x0,x1,x2,y0,y1,y2
105 real(8) c1,c2,t0,t1,t2
106 ! evaluate the polynomial coefficients
107 x0=xa(1)
108 x1=xa(2)-x0; x2=xa(3)-x0
109 y0=ya(1)
110 y1=ya(2)-y0; y2=ya(3)-y0
111 t0=1.d0/(x1*x2*(x2-x1))
112 t1=x1*y2; t2=x2*y1
113 c1=x2*t2-x1*t1
114 c2=t1-t2
115 t1=x-x0
116 ! evaluate the polynomial
117 poly3=y0+t0*t1*(c1+c2*t1)
118 end function
119 
120 pure real(8) function poly4i(xa,ya,x)
121 implicit none
122 ! arguments
123 real(8), intent(in) :: xa(4),ya(4),x
124 ! local variables
125 real(8) x0,x1,x2,x3,y0,y1,y2,y3
126 real(8) c1,c2,c3,t0,t1,t2,t3,t4,t5,t6
127 ! evaluate the polynomial coefficients
128 x0=xa(1)
129 x1=xa(2)-x0; x2=xa(3)-x0; x3=xa(4)-x0
130 y0=ya(1)
131 y1=ya(2)-y0; y2=ya(3)-y0; y3=ya(4)-y0
132 t4=x1-x2; t5=x1-x3; t6=x2-x3
133 t1=x1*x2*y3; t2=x2*x3*y1; t3=x1*x3
134 t0=1.d0/(x2*t3*t4*t5*t6)
135 t3=t3*y2
136 c3=t1*t4+t2*t6-t3*t5
137 t4=x1**2; t5=x2**2; t6=x3**2
138 c2=t1*(t5-t4)+t2*(t6-t5)+t3*(t4-t6)
139 c1=t1*(x2*t4-x1*t5)+t2*(x3*t5-x2*t6)+t3*(x1*t6-x3*t4)
140 t1=x-x0
141 ! integrate the polynomial
142 poly4i=t1*(y0+t0*t1*(0.5d0*c1+t1*(0.3333333333333333333d0*c2+0.25d0*c3*t1)))
143 end function
144 
145 end subroutine
146 !EOC
147 
pure real(8) function poly3(xa, ya, x)
Definition: rdiracint.f90:101
pure subroutine rschrodint(sol, l, e, nr, r, vr, nn, p0, p1, q0, q1)
Definition: rschrodint.f90:10
real(8), parameter sol
Definition: modmain.f90:1251
pure real(8) function poly4i(xa, ya, x)
Definition: rdiracint.f90:122