The Elk Code
 
Loading...
Searching...
No Matches
rdirac.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 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: rdirac
8! !INTERFACE:
9subroutine rdirac(sol,n,l,k,nr,r,vr,eval,g0,f0)
10! !INPUT/OUTPUT PARAMETERS:
11! sol : speed of light in atomic units (in,real)
12! n : principal quantum number (in,integer)
13! l : quantum number l (in,integer)
14! k : quantum number k (l or l+1) (in,integer)
15! nr : number of radial mesh points (in,integer)
16! r : radial mesh (in,real(nr))
17! vr : potential on radial mesh (in,real(nr))
18! eval : eigenvalue without rest-mass energy (inout,real)
19! g0 : major component of the radial wavefunction (out,real(nr))
20! f0 : minor component of the radial wavefunction (out,real(nr))
21! !DESCRIPTION:
22! Finds the solution to the radial Dirac equation for a given potential $v(r)$
23! and quantum numbers $n$, $k$ and $l$. The method involves integrating the
24! equation using the predictor-corrector method and adjusting $E$ until the
25! number of nodes in the wavefunction equals $n-l-1$. The calling routine must
26! provide an initial estimate for the eigenvalue. Note that the arrays
27! {\tt g0} and {\tt f0} represent the radial functions multiplied by $r$.
28!
29! !REVISION HISTORY:
30! Created September 2002 (JKD)
31!EOP
32!BOC
33implicit none
34! arguments
35real(8), intent(in) :: sol
36integer, intent(in) :: n,l,k,nr
37real(8), intent(in) :: r(nr),vr(nr)
38real(8), intent(inout) :: eval
39real(8), intent(out) :: g0(nr),f0(nr)
40! local variables
41integer, parameter :: maxit=2000
42integer kpa,it,ir
43integer nn,nnd,nndp
44! energy convergence tolerance
45real(8), parameter :: eps=1.d-12
46real(8) t1,de
47! automatic arrays
48real(8) g1(nr),f1(nr),fr(nr)
49! external functions
50real(8), external :: splint
51if (k < 1) then
52 write(*,*)
53 write(*,'("Error(rdirac): k < 1 : ",I8)') k
54 write(*,*)
55 stop
56end if
57if (k > n) then
58 write(*,*)
59 write(*,'("Error(rdirac): incompatible n and k : ",2I8)') n,k
60 write(*,*)
61 stop
62end if
63if ((k == n).and.(l /= k-1)) then
64 write(*,*)
65 write(*,'("Error(rdirac): incompatible n, k and l : ",3I8)') n,k,l
66 write(*,*)
67 stop
68end if
69if (k == l) then
70 kpa=k
71else if (k == l+1) then
72 kpa=-k
73else
74 write(*,*)
75 write(*,'("Error(rdirac): incompatible l and k : ",2I8)') l,k
76 write(*,*)
77 stop
78end if
79if (nr < 4) then
80 write(*,*)
81 write(*,'("Error(rdirac): nr < 4 : ",I8)') nr
82 write(*,*)
83 stop
84end if
85de=1.d0
86nndp=0
87do it=1,maxit
88! integrate the Dirac equation
89 call rdiracint(sol,kpa,eval,nr,r,vr,nn,g0,g1,f0,f1)
90! check the number of nodes
91 nnd=nn-(n-l-1)
92 if (nnd > 0) then
93 eval=eval-de
94 else
95 eval=eval+de
96 end if
97 if (it > 1) then
98 if ((nnd /= 0).or.(nndp /= 0)) then
99 if (nnd*nndp < 1) then
100 de=de*0.5d0
101 else
102 de=de*1.1d0
103 end if
104 end if
105 end if
106 nndp=nnd
107 if (de < eps*(abs(eval)+1.d0)) goto 10
108end do
109write(*,*)
110write(*,'("Warning(rdirac): maximum iterations exceeded")')
11110 continue
112! find effective infinity and set wavefunction to zero after that point
113! major component
114do ir=nr,2,-1
115 if ((g0(ir-1)*g0(ir) < 0.d0).or.(g1(ir-1)*g1(ir) < 0.d0)) then
116 g0(ir:nr)=0.d0
117 exit
118 end if
119end do
120! minor component
121do ir=nr,2,-1
122 if ((f0(ir-1)*f0(ir) < 0.d0).or.(f1(ir-1)*f1(ir) < 0.d0)) then
123 f0(ir:nr)=0.d0
124 exit
125 end if
126end do
127! normalise
128do ir=1,nr
129 fr(ir)=g0(ir)**2+f0(ir)**2
130end do
131t1=splint(nr,r,fr)
132t1=sqrt(abs(t1))
133if (t1 <= 0.d0) then
134 write(*,*)
135 write(*,'("Error(rdirac): zero wavefunction")')
136 write(*,*)
137 stop
138end if
139t1=1.d0/t1
140g0(1:nr)=t1*g0(1:nr)
141f0(1:nr)=t1*f0(1:nr)
142end subroutine
143!EOC
144
subroutine rdirac(sol, n, l, k, nr, r, vr, eval, g0, f0)
Definition rdirac.f90:10
pure subroutine rdiracint(sol, kpa, e, nr, r, vr, nn, g0, g1, f0, f1)
Definition rdiracint.f90:10
real(8) function splint(n, x, f)
Definition splint.f90:7