25 use,
intrinsic :: iso_fortran_env
51 integer :: no_l_channels
52 real(real64),
allocatable :: vps(:, :)
53 real(real64),
allocatable :: KB(:,:)
54 real(real64),
allocatable :: dkbcos(:)
55 real(real64),
allocatable,
private :: dknorm(:)
56 real(real64),
allocatable :: kb_radius(:)
57 integer,
allocatable :: kb_nrval(:)
59 integer :: so_no_l_channels
60 real(real64),
allocatable :: so_vps(:,:)
61 real(real64),
allocatable,
private :: so_KB(:,:)
62 real(real64),
allocatable,
private :: so_dkbcos(:)
63 real(real64),
allocatable,
private :: so_dknorm(:)
64 real(real64),
allocatable,
private :: so_kb_radius(:)
65 integer,
allocatable :: so_kb_nrval(:)
67 real(real64),
allocatable :: vlocal(:)
68 real(real64),
allocatable :: rphi(:, :,:)
70 logical :: core_corrections
71 real(real64),
allocatable :: chcore(:)
78 type(ps_in_grid_t),
intent(out) :: ps
79 integer,
intent(in) :: flavor, nrval
80 real(real64),
intent(in) :: a, b
81 integer,
intent(in) :: no_l, so_no_l
89 ps% no_l_channels = no_l
90 ps%so_no_l_channels = so_no_l
93 safe_allocate(ps%vps (1:ps%g%nrval, 1:no_l))
94 safe_allocate(ps%chcore (1:ps%g%nrval))
95 safe_allocate(ps%vlocal (1:ps%g%nrval))
97 safe_allocate(ps%rphi (1:ps%g%nrval, 1:no_l, 1:3))
98 safe_allocate(ps%KB (1:ps%g%nrval, 1:no_l))
99 safe_allocate(ps%dkbcos(1:no_l))
100 safe_allocate(ps%dknorm(1:no_l))
101 safe_allocate(ps%kb_radius(1:no_l+1))
102 safe_allocate(ps%kb_nrval(1:no_l+1))
104 if (so_no_l > 0)
then
105 safe_allocate(ps%so_vps(1:ps%g%nrval, 1:so_no_l))
106 safe_allocate(ps%so_KB (1:ps%g%nrval, 1:so_no_l))
107 safe_allocate(ps%so_dkbcos(1:so_no_l))
108 safe_allocate(ps%so_dknorm(1:so_no_l))
109 safe_allocate(ps%so_kb_radius(1:so_no_l))
110 safe_allocate(ps%so_kb_nrval(1:so_no_l))
119 type(ps_in_grid_t),
intent(inout) :: ps
123 safe_deallocate_a(ps%vps)
124 safe_deallocate_a(ps%chcore)
125 safe_deallocate_a(ps%vlocal)
127 safe_deallocate_a(ps%rphi)
128 safe_deallocate_a(ps%KB)
129 safe_deallocate_a(ps%dkbcos)
130 safe_deallocate_a(ps%dknorm)
131 safe_deallocate_a(ps%kb_radius)
132 safe_deallocate_a(ps%kb_nrval)
134 if (ps%so_no_l_channels > 0)
then
135 safe_deallocate_a(ps%so_vps)
136 safe_deallocate_a(ps%so_KB)
137 safe_deallocate_a(ps%so_dkbcos)
138 safe_deallocate_a(ps%so_dknorm)
139 safe_deallocate_a(ps%so_kb_radius)
140 safe_deallocate_a(ps%so_kb_nrval)
152 integer,
intent(in) :: l_loc
153 real(real64),
intent(in) :: rcore
157 real(real64) :: a, b, qtot
158 real(real64),
allocatable :: rho(:)
163 ps%vlocal(:) = ps%vps(:, l_loc+1)
165 else if (l_loc == -1)
then
167 message(1) =
"For the moment, Vanderbilt local potentials are only possible with tm grids."
171 a = 1.82_real64 / rcore
173 safe_allocate(rho(1:ps%g%nrval))
175 do ir = 1, ps%g%nrval
176 rho(ir) =
exp(-(
sinh(a*b*ps%g%rofi(ir)) /
sinh(b))**2)
177 rho(ir) =
m_four *
m_pi * rho(ir) * ps%g%r2ofi(ir)
180 qtot = sum(rho(:)*ps%g%drdi(:))
181 rho(:) = - rho(:)*(ps%zval/qtot)
183 call vhrtre(rho, ps%vlocal, ps%g%rofi, ps%g%drdi, ps%g%s, ps%g%nrval, ps%g%a)
184 ps%vlocal(1) = ps%vlocal(2)
186 safe_deallocate_a(rho)
203 do l = 1, ps%no_l_channels
204 ps%KB(2:, l) = (ps%vps(2:, l) - ps%vlocal(2:))*(ps%rphi(2:, l, 1)/ps%g%rofi(2:))*ps%dknorm(l)
209 do l = 1, ps%so_no_l_channels
210 ps%so_KB(2:, l) = ps%so_vps(2:, l)*(ps%rphi(2:, l, 1)/ps%g%rofi(2:))*ps%dknorm(l)
227 integer,
intent(in) :: lloc
230 real(real64) :: dnrm, avgv, vphi
231 real(real64),
parameter :: tol = 1.0e-20_real64
235 do l = 1, ps%no_l_channels
236 if (l-1 == lloc)
then
244 do ir = 1, ps%g%nrval
245 vphi = (ps%vps(ir, l) - ps%vlocal(ir))*ps%rphi(ir, l, 1)
246 dnrm = dnrm + vphi*vphi*ps%g%drdi(ir)
247 avgv = avgv + vphi*ps%rphi(ir, l, 1)*ps%g%drdi(ir)
249 ps%dkbcos(l) = dnrm/safe_tol(avgv, tol)
250 ps%dknorm(l) =
m_one/safe_tol(
sqrt(dnrm), tol)
253 do l = 1, ps%so_no_l_channels
256 do ir = 1, ps%g%nrval
257 vphi = ps%so_vps(ir, l)*ps%rphi(ir, l, 1)
258 dnrm = dnrm + vphi*vphi*ps%g%drdi(ir)
259 avgv = avgv + vphi*ps%rphi(ir, l, 1)*ps%g%drdi(ir)
261 ps%so_dkbcos(l) = dnrm/safe_tol(avgv, tol)
262 ps%so_dknorm(l) =
m_one/safe_tol(
sqrt(dnrm), tol)
272 integer,
intent(in) :: lloc
275 real(real64) :: dincv, phi
276 real(real64),
parameter :: threshold = 1.0e-6_real64
281 do ir = ps%g%nrval-1, 2, -1
282 dincv = abs(ps%vlocal(ir)*ps%g%rofi(ir) +
m_two*ps%zval)
283 if (dincv > threshold)
exit
285 ps%kb_radius(ps%no_l_channels+1) = ps%g%rofi(ir + 1)
286 ps%kb_nrval(ps%no_l_channels+1) = ir + 1
289 do l = 1, ps%no_l_channels
290 if (l-1 == lloc)
then
296 do ir = ps%g%nrval-1, 2, -1
297 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%dknorm(l)
298 dincv = abs((ps%vps(ir, l) - ps%vlocal(ir))*phi)
300 if (dincv > threshold)
exit
302 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%dknorm(l)
305 ps%kb_radius(l) = ps%g%rofi(ir + 1)
306 ps%kb_nrval(l) = ir + 1
310 do l = 1, ps%so_no_l_channels
311 do ir = ps%g%nrval-1, 2, -1
312 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%so_dknorm(l)
313 dincv = abs((ps%so_vps(ir, l))*phi)
315 if (dincv > threshold)
exit
318 ps%so_kb_radius(l) = ps%g%rofi(ir + 1)
319 ps%so_kb_nrval(l) = ir + 1
338 do l = 1, ps%no_l_channels
339 nrm =
sqrt(sum(ps%g%drdi(:)*ps%rphi(:, l, 1)**2))
340 nrm = abs(nrm -
m_one)
341 if (nrm > 1.0e-5_real64)
then
342 write(
message(1),
'(a,i2,a)')
"Eigenstate for l = ", l-1,
' is not normalized.'
343 write(
message(2),
'(a, f12.6,a)')
'(abs(1 - norm) = ', nrm,
')'
354 real(real64),
intent(in) :: x(:)
355 real(real64),
intent(in) :: y(:)
356 logical,
optional,
intent(in) :: high_order
358 real(real64) :: x1, x2, x3
359 real(real64) :: y1, y2, y3
370 y0 = y1*x2*x3*(x2 - x3) + y2*x1*x3*(x3 - x1) + y3*x1*x2*(x1 - x2)
371 y0 = y0/((x1 - x2)*(x1 - x3)*(x2 - x3))
375 y0 = y1 - (y2 - y1)*x1/(x2 - x1)
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double sinh(double __x) __attribute__((__nothrow__
subroutine, public vhrtre(rho, v, r, drdi, srdrdi, nr, a)
VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY....
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_one
integer, parameter, public logrid_psf
log grid used in Troullier-Martins code
subroutine, public logrid_end(grid)
subroutine, public logrid_init(grid, flavor, aa, bb, nrval)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public ps_in_grid_kb_projectors(ps)
KB-projectors kb = (vps - vlocal) |phi> * dknorm.
subroutine, public ps_in_grid_end(ps)
subroutine, public ps_in_grid_cutoff_radii(ps, lloc)
subroutine, public ps_in_grid_init(ps, flavor, a, b, nrval, no_l, so_no_l)
subroutine, public ps_in_grid_kb_cosines(ps, lloc)
KB-cosines and KB-norms: dkbcos stores the KB "cosines:" || (v_l - v_local) phi_l ||^2 / < (v_l - v_l...
subroutine, public ps_in_grid_vlocal(ps, l_loc, rcore, namespace)
real(real64) function, public first_point_extrapolate(x, y, high_order)
subroutine, public ps_in_grid_check_rphi(ps, namespace)
checks normalization of the pseudo wavefunctions