24 use,
intrinsic :: iso_fortran_env
42 integer,
parameter,
public :: &
43 LOGRID_PSF = 1, & !< log grid used in Troullier-Martins code
53 real(real64),
allocatable :: rofi(:)
54 real(real64),
allocatable :: r2ofi(:)
55 real(real64),
allocatable :: drdi(:)
56 real(real64),
allocatable :: s(:)
63 type(logrid_t),
intent(out) :: grid
64 integer,
intent(in) :: flavor
65 real(real64),
intent(in) :: aa, bb
66 integer,
intent(in) :: nrval
68 real(real64) :: rpb, ea
73 assert(flavor == logrid_psf .or. flavor ==
logrid_cpi)
83 safe_allocate(grid%rofi(1:nrval))
84 safe_allocate(grid%r2ofi(1:nrval))
85 safe_allocate(grid%drdi(1:nrval))
86 safe_allocate(grid%s(1:nrval))
88 select case (grid%flavor)
93 grid%drdi(ir) = aa*rpb
95 grid%rofi(ir) = bb*(
exp(aa*(ir-1)) -
m_one)
104 grid%drdi(2) = bb*rpb
105 do ir = 3, grid%nrval
106 grid%rofi(ir) = grid%rofi(ir-1)*aa
107 grid%drdi(ir) = grid%rofi(ir)*rpb
112 do ir = 1, grid%nrval
113 grid%s(ir) =
sqrt(grid%drdi(ir))
114 grid%r2ofi(ir) = grid%rofi(ir)**2
123 type(logrid_t),
intent(inout) :: grid
127 safe_deallocate_a(grid%rofi)
128 safe_deallocate_a(grid%r2ofi)
129 safe_deallocate_a(grid%drdi)
130 safe_deallocate_a(grid%s)
138 type(
logrid_t),
intent(in) :: grid_in
139 type(logrid_t),
intent(inout) :: grid_out
145 grid_out%flavor = grid_in%flavor
146 grid_out%a = grid_in%a
147 grid_out%b = grid_in%b
148 grid_out%nrval = grid_in%nrval
150 safe_allocate(grid_out%rofi (1:grid_out%nrval))
151 safe_allocate(grid_out%r2ofi(1:grid_out%nrval))
152 safe_allocate(grid_out%drdi (1:grid_out%nrval))
153 safe_allocate(grid_out%s (1:grid_out%nrval))
155 grid_out%rofi(:) = grid_in%rofi(:)
156 grid_out%r2ofi(:) = grid_in%r2ofi(:)
157 grid_out%drdi(:) = grid_in%drdi(:)
158 grid_out%s(:) = grid_in%s(:)
167 real(real64),
intent(in) :: rofi
174 do ir = 1, grid%nrval-1
176 if (rofi >= grid%rofi(ir).and.rofi < grid%rofi(ir+1))
then
177 if (abs(rofi-grid%rofi(ir)) < abs(rofi-grid%rofi(ir+1)))
then
194 real(real64),
intent(in) :: ff(:)
195 real(real64),
intent(out) :: dfdr(:)
201 dfdr(1) = (ff(2) - ff(1))/(grid%rofi(2) - grid%rofi(1))
202 do ii = 2, grid%nrval-1
203 dfdr(ii) = (ff(ii+1) - ff(ii-1))/(grid%rofi(ii+1) - grid%rofi(ii-1))
205 dfdr(grid%nrval) = (ff(grid%nrval) - ff(grid%nrval-1))/(grid%rofi(grid%nrval) - grid%rofi(grid%nrval-1))
211 real(real64)
pure function
logrid_radius(grid) result(radius)
214 radius = grid%rofi(grid%nrval)
219 type(namespace_t),
intent(in) :: namespace
220 integer,
intent(in) :: zz
221 real(real64),
intent(out) :: aa, bb
222 integer,
intent(out) :: np
224 real(real64) :: xmin, xmax, a1, a2, f1, fm
233 xmin =
sqrt(real(zz, real64) )*1e-5_real64
234 xmax =
sqrt(real(zz, real64) )*30.0_real64
235 np =
floor(
sqrt(real(zz, real64) )*200_real64)
237 np =
floor(np/m_two)*2+1
240 f1 =
func(xmin, xmax, real(np, real64) , a1)
243 aa = (a2 + a1)*m_half
244 fm =
func(xmin, xmax, real(np, real64) , aa)
245 if (m_half*abs(a1 - a2) < 1.0e-16_real64)
exit
246 if (fm*f1 > m_zero)
then
254 bb = xmin/(
exp(aa)-m_one)
256 write(message(1),
'(a,es13.6,a,es13.6,a,i4)')
'Debug: Log grid parameters: a = ', aa, &
257 ' b = ', bb,
' np = ', np
258 call messages_info(1, namespace=namespace, debug_only=.
true.)
262 real(real64) function
func(r1, rn, n, a)
263 real(real64),
intent(in) :: r1, rn, a, n
264 if((n-m_one)*a < m_max_exp_arg)
then
265 func =
exp((n-m_one)*a)*r1 - m_one*r1 - rn*
exp(a) + rn*m_one
267 func = m_huge*r1 - m_one*r1 - rn*
exp(a) + rn*m_one
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
real(real64) function func(r1, rn, n, a)
real(real64), parameter, public m_zero
real(real64), parameter, public m_one
subroutine, public logrid_find_parameters(namespace, zz, aa, bb, np)
real(real64) pure function, public logrid_radius(grid)
subroutine, public logrid_copy(grid_in, grid_out)
integer, parameter, public logrid_cpi
log grid used in FHI code
integer function, public logrid_index(grid, rofi)
subroutine, public logrid_end(grid)
subroutine, public derivative_in_log_grid(grid, ff, dfdr)
subroutine, public logrid_init(grid, flavor, aa, bb, nrval)