The Elk Code
findband.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: findband
8 ! !INTERFACE:
9 subroutine findband(sol,l,nr,r,vr,eps,demax,e,fnd)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! sol : speed of light in atomic units (in,real)
12 ! l : angular momentum quantum number (in,integer)
13 ! nr : number of radial mesh points (in,integer)
14 ! r : radial mesh (in,real(nr))
15 ! vr : potential on radial mesh (in,real(nr))
16 ! eps : energy search tolerance (in,real)
17 ! demax : maximum allowed change from the input energy; enforced only if e < 0
18 ! (in,real)
19 ! e : input energy and returned band energy (inout,real)
20 ! fnd : set to .true. if the band energy is found (out,logical)
21 ! !DESCRIPTION:
22 ! Finds the band energies for a given radial potential and angular momentum.
23 ! This is done by first searching upwards in energy, starting from the input
24 ! energy plus the offset energy, until the radial wavefunction at the
25 ! muffin-tin radius is zero. This is the energy at the top of the band,
26 ! denoted $E_{\rm t}$. A downward search is now performed from $E_{\rm t}$
27 ! until the slope of the radial wavefunction at the muffin-tin radius is zero.
28 ! This energy, $E_{\rm b}$, is at the bottom of the band. The band energy is
29 ! taken as $(E_{\rm t}+E_{\rm b})/2$. If either $E_{\rm t}$ or $E_{\rm b}$
30 ! cannot be found then the band energy is set to the input value.
31 !
32 ! !REVISION HISTORY:
33 ! Created September 2004 (JKD)
34 ! Added two-pass loop, October 2013 (JKD)
35 !EOP
36 !BOC
37 implicit none
38 ! arguments
39 real(8), intent(in) :: sol
40 integer, intent(in) :: l,nr
41 real(8), intent(in) :: r(nr),vr(nr),eps,demax
42 real(8), intent(inout) :: e
43 logical, intent(out) :: fnd
44 ! local variables
45 logical ft,fb
46 ! maximum number of steps
47 integer, parameter :: maxstp=250
48 integer ip,ie,nn
49 ! initial step size
50 real(8), parameter :: de0=0.001d0
51 real(8) de,et,eb,t,tp
52 ! automatic arrays
53 real(8) p0(nr),p1(nr),q0(nr),q1(nr)
54 ft=.false.
55 fb=.false.
56 fnd=.false.
57 et=e
58 eb=e
59 ! two-pass loop
60 do ip=1,2
61 ! find the top of the band
62  tp=0.d0
63  de=de0
64  do ie=1,maxstp
65  et=et+de
66  if (e < 0.d0) then
67  if (et > e+demax) exit
68  end if
69  call rschrodint(sol,l,et,nr,r,vr,nn,p0,p1,q0,q1)
70  t=p0(nr)
71  if (ie > 1) then
72  if (t*tp <= 0.d0) then
73  if (abs(de) < eps) then
74  if (fb) goto 10
75  ft=.true.
76  eb=et+5.d0*abs(de0)
77  exit
78  end if
79  de=-0.5d0*de
80  else
81  de=1.5d0*de
82  end if
83  end if
84  tp=t
85  end do
86  if (fb) return
87 ! find the bottom of the band
88  tp=0.d0
89  de=-de0
90  do ie=1,maxstp
91  eb=eb+de
92  if (eb < e-demax) return
93  call rschrodint(sol,l,eb,nr,r,vr,nn,p0,p1,q0,q1)
94  t=p1(nr)
95  if (ie > 1) then
96  if (t*tp <= 0.d0) then
97  if (abs(de) < eps) then
98  if (ft) goto 10
99  fb=.true.
100  et=eb-5.d0*abs(de0)
101  exit
102  end if
103  de=-0.5d0*de
104  else
105  de=1.5d0*de
106  end if
107  end if
108  tp=t
109  end do
110 end do
111 return
112 10 continue
113 ! set the band energy halfway between top and bottom
114 e=(et+eb)/2.d0
115 fnd=.true.
116 end subroutine
117 !EOC
118 
pure subroutine rschrodint(sol, l, e, nr, r, vr, nn, p0, p1, q0, q1)
Definition: rschrodint.f90:10
subroutine findband(sol, l, nr, r, vr, eps, demax, e, fnd)
Definition: findband.f90:10