The Elk Code
findlambda.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2009 F. Bultmark, F. Cricchio and L. Nordstrom.
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: findlambda
8 ! !INTERFACE:
9 subroutine findlambda(is,l,ufix,lambda0,lambda)
10 use modmpi
11 ! !INPUT/OUTPUT PARAMETERS:
12 ! is : species type (in,integer)
13 ! l : angular momentum (in,integer)
14 ! ufix : fixed U (in,integer)
15 ! lambda0 : starting value for screening length (inout,real)
16 ! lambda : screening length corresponding to fixed U (out,real)
17 ! !DESCRIPTION:
18 ! Find the screening length corresponding to a fixed value of $U$ by using the
19 ! half-interval method in the first few steps and then the more efficient
20 ! secant method. For $U=0$ the code automatically sets the screening length to
21 ! ${\rm lambdamax}=50$. This value is enough to get $F^{(k)}\sim 10^{-3}$
22 ! corresponding to $U\sim 0$ (that perfectly mimics a bare DFT calculation).
23 !
24 ! !REVISION HISTORY:
25 ! Created July 2009 (Francesco Cricchio)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: is,l
31 real(8), intent(in) :: ufix
32 real(8), intent(inout) :: lambda0
33 real(8), intent(out) :: lambda
34 ! local variables
35 ! max iterations in secant algorithm
36 integer, parameter :: maxit=100
37 integer it,nit
38 ! if ufix < umin, assume lambda=lambdamax
39 real(8), parameter :: umin=1.d-4
40 ! if lambda < lambdamin, perform unscreened calculation
41 real(8), parameter :: lambdamin=1.d-2
42 ! max value of lambda
43 ! lambdamax=50 is enough to get F^(k)~1.d-3 corresponding to U~0
44 real(8), parameter :: lambdamax=50.d0
45 ! initial step for lambda
46 real(8), parameter :: dl0=0.5d0
47 real(8) f,fp,lambdap,dl,tol
48 ! external functions
49 real(8), external :: fyukawa,fyukawa0
50 ! small U limit
51 if (ufix < umin) then
52  lambda=lambdamax
53  if (mp_mpi) write(*,'("Info(findlambda): lambda set to lambdamax : ",&
54  &G18.10)') lambdamax
55  return
56 end if
57 ! first perform a search of lambda with half-interval method and low accuracy
58 ! initialize values and search upward from lambda0
59 lambda=lambda0
60 dl=dl0
61 fp=1.d0
62 tol=1.d-1
63 nit=0
64 do it=1,maxit
65  if (lambda < lambdamin) then
66 ! unscreened Slater parameters
67  f=fyukawa0(is,l,0)-ufix
68  else
69 ! screened Slater parameters
70  f=fyukawa(is,l,0,lambda)-ufix
71  end if
72  if ((f*fp) < 0) dl=-0.5d0*dl
73  lambdap=lambda
74  lambda=lambda+dl
75  fp=f
76  nit=nit+1
77  if (abs(f) < tol) goto 10
78 end do
79 10 continue
80 ! use the found value of lambda to continue the search with secant algorithm and
81 ! higher accuracy
82 tol=1.d-8
83 ! calculate F^(0)-ufix at lambdap value
84 if (lambdap < lambdamin) then
85 ! unscreened Slater parameters
86  fp=fyukawa0(is,l,0)-ufix
87 else
88 ! screened Slater parameters
89  fp=fyukawa(is,l,0,lambdap)-ufix
90 end if
91 ! start secant algorithm
92 do it=1,maxit
93 ! calculate F^(0)-ufix
94  if (lambda < lambdamin) then
95 ! unscreened Slater parameters
96  f=fyukawa0(is,l,0)-ufix
97  else
98 ! screened Slater parameters
99  f=fyukawa(is,l,0,lambda)-ufix
100  end if
101 ! if requested tolerance has been reached exit the loop
102  if (abs(f) < tol) goto 20
103 ! update lambda with secant algorithm and roll values
104  dl=-f*((lambda-lambdap)/(f-fp))
105  lambdap=lambda
106  lambda=lambda+dl
107  fp=f
108  nit=nit+1
109 end do
110 20 continue
111 if (nit >= maxit) then
112  write(*,*)
113  write(*,'("Error(findlambda): max number of iterations to obtain lambda &
114  &reached")')
115  write(*,*)
116  stop
117 else
118 ! update initial value for lambda for the next iteration in the SC loop
119 ! 0.5*dl0 is enough
120  lambda0=lambda-0.5d0*dl0
121  if (mp_mpi) write(*,'("Info(findlambda): lambda obtained in ",I4,&
122  &" iterations")') nit
123 end if
124 end subroutine
125 !EOC
126 
logical mp_mpi
Definition: modmpi.f90:17
subroutine findlambda(is, l, ufix, lambda0, lambda)
Definition: findlambda.f90:10
Definition: modmpi.f90:6