30 use,
intrinsic :: iso_fortran_env
54 character(len=256),
private :: title
55 integer,
private :: pspdat
56 integer,
private :: pspcod
57 integer,
private :: pspxc
59 integer,
private :: mmax
60 real(real64),
private :: r2well
63 character(len=5),
private :: atom_name
65 real(real64) :: rlocal
66 real(real64),
private :: rc(0:3)
67 real(real64),
private :: c(1:4)
68 real(real64) :: h(0:3, 1:3, 1:3)
69 real(real64) :: k(0:3, 1:3, 1:3)
71 type(valconf_t) :: conf
74 real(real64),
allocatable :: vlocal(:)
75 real(real64),
allocatable :: kb(:,:,:)
76 real(real64),
allocatable :: kbr(:)
77 real(real64),
allocatable :: rphi(:, :)
78 real(real64),
allocatable,
private :: eigen(:)
84 real(real64),
parameter :: eps = 1.0e-8_real64
89 subroutine hgh_init(psp, filename, namespace)
90 type(ps_hgh_t),
intent(inout) :: psp
91 character(len=*),
intent(in) :: filename
92 type(namespace_t),
intent(in) :: namespace
94 integer :: iunit, i, np
95 real(real64) :: aa, bb
99 iunit =
io_open(trim(filename), action=
'read', form=
'formatted', status=
'old')
109 do while(psp%rc(psp%l_max) > 0.01_real64)
110 psp%l_max = psp%l_max + 1
111 if (psp%l_max > 3)
exit
113 psp%l_max = psp%l_max - 1
120 safe_allocate(psp%vlocal(1:psp%g%nrval))
122 if (psp%l_max >= 0)
then
123 safe_allocate(psp%kbr(0:psp%l_max))
124 safe_allocate(psp%kb(1:psp%g%nrval, 0:psp%l_max, 1:3))
135 type(ps_hgh_t),
intent(inout) :: psp
139 if (psp%l_max >= 0)
then
140 safe_deallocate_a(psp%kbr)
141 safe_deallocate_a(psp%kb)
143 safe_deallocate_a(psp%vlocal)
144 safe_deallocate_a(psp%rphi)
145 safe_deallocate_a(psp%eigen)
154 type(
ps_hgh_t),
intent(inout) :: psp
157 integer :: l, i, ierr
161 safe_allocate(psp%rphi(1:psp%g%nrval, 1:psp%conf%p))
162 safe_allocate(psp%eigen(1:psp%conf%p))
181 write(
message(1),
'(a)')
'The algorithm that calculates atomic wavefunctions could not'
182 write(
message(2),
'(a)')
'do its job. The program will continue, but expect poor'
183 write(
message(3),
'(a)')
'convergence properties.'
197 real(real64),
intent(out) :: eigen(:, :)
204 eigen(i, :) = psp%eigen(i)
212 integer,
intent(in) :: unit
213 type(
ps_hgh_t),
intent(out) :: params
219 integer :: i, iostat, j, k, nn, nnonloc, lloc
220 character(len=VALCONF_STRING_LENGTH) :: line
231 read(unit, *) params%title
232 read(unit, *) params%conf%z, params%z_val, params%pspdat
233 read(unit, *, iostat=iostat) params%pspcod, params%pspxc, params%l_max, lloc, params%mmax, params%r2well
235 select case (params%pspcod)
241 read(unit,
'(a)') line
242 do while((iostat /= 0) .and. (j > 0))
244 read(line, *, iostat=iostat) params%rlocal, params%c(1:j)
246 if (j<1)
read(line, *, iostat=iostat) params%atom_name, params%z_val, params%rlocal
247 if (iostat /= 0)
then
253 read(unit,
'(a)', iostat = iostat) line
254 if (iostat /= 0)
then
261 do while((iostat /= 0) .and. (j > 0))
263 read(line, *, iostat=iostat) params%rc(0), (params%h(0, i, i), i = 1, j)
272 read(unit,
'(a)', iostat = iostat) line
273 if (iostat /= 0)
exit kloop
276 do while((iostat /= 0) .and. (j > 0))
278 read(line, *, iostat = iostat) params%rc(k), (params%h(k, i, i), i = 1, j)
280 if (abs(params%rc(k)) <=
m_epsilon)
exit kloop
281 read(unit,
'(a)') line
284 do while((iostat /= 0) .and. (j>0))
286 read(line, *, iostat = iostat) (params%k(k, i, i), i = 1, 3)
293 params%h(0, 2, 3) = -
m_half *
sqrt(100.0_real64/63.0_real64) * params%h(0, 3, 3)
295 params%h(1, 1, 3) =
m_one/6.0_real64 *
sqrt(35.0_real64/11.0_real64) * params%h(1, 3, 3)
296 params%h(1, 2, 3) = -
m_one/6.0_real64 * (14.0_real64 /
sqrt(11.0_real64)) * params%h(1, 3, 3)
297 params%h(2, 1, 2) = -
m_half *
sqrt(7.0_real64/9.0_real64) * params%h(2, 2, 2)
298 params%h(2, 1, 3) =
m_half *
sqrt(63.0_real64/143.0_real64) * params%h(2, 3, 3)
299 params%h(2, 2, 3) = -
m_half * (18.0_real64/
sqrt(143.0_real64)) * params%h(2, 3, 3)
303 params%k(0, 2, 3) = -
m_half *
sqrt(100.0_real64/63.0_real64) * params%k(0, 3, 3)
305 params%k(1, 1, 3) =
m_one/6.0_real64 *
sqrt(35.0_real64/11.0_real64) * params%k(1, 3, 3)
306 params%k(1, 2, 3) = -
m_one/6.0_real64 * (14.0_real64 /
sqrt(11.0_real64)) * params%k(1, 3, 3)
307 params%k(2, 1, 2) = -
m_half *
sqrt(7.0_real64/9.0_real64) * params%k(2, 2, 2)
308 params%k(2, 1, 3) =
m_half *
sqrt(63.0_real64/143.0_real64) * params%k(2, 3, 3)
309 params%k(2, 2, 3) = -
m_half * (18.0_real64/
sqrt(143.0_real64)) * params%k(2, 3, 3)
314 read(unit, *) params%rlocal, nn, params%c(1)
315 read(unit, *) nnonloc
316 read(unit,
'(a)', iostat = iostat) line
318 if (iostat /= 0)
exit kloop2
321 do while((iostat /= 0) .and. (j > 0))
323 read(line, *, iostat = iostat) params%rc(k), nn, (params%h(k, 1, i), i = 1, j)
325 if (abs(params%rc(k)) <=
m_epsilon)
exit kloop2
326 read(unit,
'(a)') line
329 read(line, *) (params%h(1, 2, i), i = 2, 3)
330 read(unit,
'(a)') line
331 read(line, *) params%h(1, 3, 3)
333 read(unit,
'(a)', iostat=iostat) line
334 if (iostat /= 0)
exit kloop2
335 read(line, *, iostat=iostat) (params%k(1, 1, i), i = 1, 3)
336 if (iostat /= 0)
continue
337 read(unit,
'(a)') line
338 read(line, *) (params%k(1, 2, i), i = 2, 3)
339 read(unit,
'(a)') line
340 read(line, *) params%k(1, 3, 3)
341 read(unit,
'(a)', iostat = iostat) line
343 read(line, *) params%h(0, 2, 2)
345 read(unit,
'(a)', iostat=iostat) line
346 if (iostat /= 0)
exit kloop2
347 read(line, *, iostat=iostat) (params%k(1, 1, i), i = 1, 2)
348 if (iostat /= 0)
continue
349 read(unit,
'(a)') line
350 read(line, *) params%k(1, 2, 2)
351 read(unit,
'(a)', iostat = iostat) line
353 message(1) =
"Error parsing the pseudopotential"
359 message(1) =
"Inconsistency in pseudopotential file:"
360 write(
message(2),
'(a,i2)')
" expecting pspcod = 3, but found ", params%pspcod
369 params%h(k, j, i) = params%h(k, i, j)
370 params%k(k, j, i) = params%k(k, i, j)
382 type(
ps_hgh_t),
intent(inout) :: psp
385 real(real64) :: dincv, tmp
386 real(real64),
parameter :: threshold = 1.0e-4_real64
393 do ir = psp%g%nrval, 2, -1
394 dincv = abs(psp%kb(ir, l, i))
395 if (dincv > threshold)
exit
397 tmp = psp%g%rofi(ir + 1)
398 psp%kbr(l) = max(tmp, psp%kbr(l))
411 real(real64),
intent(in) :: r(:)
412 integer,
intent(in) :: np
413 real(real64),
intent(inout) :: vloc(:)
416 real(real64) :: r1, r2
421 real(real64),
parameter :: tol = 5e-5_real64
435 vloc(ip) = vloc(ip) +
exp(-
m_half*r2) * (p%c(1) + r2*(p%c(2) + r2*(p%c(3) + p%c(4)*r2)))
448 real(real64),
intent(in) :: r(:)
449 integer,
intent(in) :: np
450 integer,
intent(in) :: i
451 integer,
intent(in) :: l
452 real(real64),
intent(inout) :: proj(:)
455 real(real64) :: x, y, rr, arg
459 x = l + real(4*i-1, real64) /
m_two
466 if (l /= 0 .or. i /= 1)
then
467 rr = r(ip) ** (l + 2*(i-1))
470 arg = -r(ip)**2/(
m_two*p%rc(l)**2)
473 (p%rc(l)**(l + real(4*i-1, real64) /
m_two) * x)
485 type(
ps_hgh_t),
intent(inout) :: psp
486 integer,
intent(out) :: ierr
489 integer :: iter, ir, l, nnode, nprin, i, j, irr, n, k
490 real(real64) :: vtot, a2b4, diff, nonl
491 real(real64),
allocatable :: prev(:, :), rho(:, :), ve(:, :)
492 real(real64),
parameter :: tol = 1.0e-4_real64
493 real(real64) :: e, z, dr, rmax
494 real(real64),
allocatable :: s(:), hato(:), g(:), y(:)
496 integer,
parameter :: max_iter = 200
503 safe_allocate( s(1:psp%g%nrval))
504 safe_allocate(hato(1:psp%g%nrval))
505 safe_allocate( g(1:psp%g%nrval))
506 safe_allocate( y(1:psp%g%nrval))
507 safe_allocate(prev(1:psp%g%nrval, 1))
508 safe_allocate( rho(1:psp%g%nrval, 1))
509 safe_allocate( ve(1:psp%g%nrval, 1))
517 s(2:psp%g%nrval) = psp%g%drdi(2:psp%g%nrval)*psp%g%drdi(2:psp%g%nrval)
525 self_consistent:
do while(diff > tol)
526 prev(:, :) = rho(:, :)
528 if (iter > max_iter)
then
535 do ir = 2, psp%g%nrval
536 vtot = 2*psp%vlocal(ir) + ve(ir, 1) + real(l*(l + 1), real64)/(psp%g%rofi(ir)**2)
538 if (iter>2 .and. psp%l_max >= 0 .and. abs(psp%rphi(ir, n)) > 1.0e-7_real64)
then
541 do irr = 2, psp%g%nrval
542 if(psp%kb(irr,l,j) < 1e-100_real64 .or. psp%kb(ir, l, i) < 1e-100_real64) cycle
543 nonl = nonl + psp%h(l, i, j)*psp%kb(ir, l, i)* &
544 psp%g%drdi(irr)*psp%g%rofi(irr)*psp%rphi(irr, n)*psp%kb(irr,l,j)
548 nonl =
m_two * nonl / psp%rphi(ir, n) * psp%g%rofi(ir)
551 hato(ir) = vtot*s(ir) + a2b4
557 if (psp%conf%l(k) == psp%conf%l(n)) nnode = 2
561 e = -((psp%z_val/real(nprin, real64) )**2)
568 rmax = psp%g%rofi(psp%g%nrval)
569 call egofv(namespace, hato, s, psp%g%nrval, e, g, y, l, z, &
570 real(psp%g%a, real64),
real(psp%g%b, real64), rmax, nprin, nnode, dr, ierr)
571 if (ierr /= 0)
exit self_consistent
574 psp%rphi(2:psp%g%nrval, n) = g(2:psp%g%nrval) *
sqrt(psp%g%drdi(2:psp%g%nrval))
575 psp%rphi(1, n) = psp%rphi(2, n)
579 rho(1:psp%g%nrval, 1) = rho(1:psp%g%nrval, 1) + psp%conf%occ(n,1)*psp%rphi(1:psp%g%nrval, n)**2
581 if (iter>1) rho(:, :) =
m_half*(rho(:, :) + prev(:, :))
582 diff =
sqrt(sum(psp%g%drdi(2:psp%g%nrval)*(rho(2:psp%g%nrval, 1)-prev(2:psp%g%nrval, 1))**2))
583 call atomhxc(namespace,
'LDA', psp%g, 1, rho, ve)
585 end do self_consistent
591 e =
sqrt(sum(psp%g%drdi(2:psp%g%nrval)*psp%rphi(2:psp%g%nrval, n)**2))
593 if (e > 1.0e-5_real64)
then
594 write(
message(1),
'(a,i2,a)')
"Eigenstate for n = ", n ,
' is not normalized'
595 write(
message(2),
'(a, f12.6,a)')
'(abs(1-norm) = ', e,
')'
601 psp%eigen = psp%eigen /
m_two
606 safe_deallocate_a(hato)
609 safe_deallocate_a(prev)
610 safe_deallocate_a(rho)
611 safe_deallocate_a(ve)
618 subroutine hgh_debug(psp, dir, namespace)
620 character(len=*),
intent(in) :: dir
623 integer :: hgh_unit, loc_unit, dat_unit, kbp_unit, wav_unit, i, l, k
624 character(len=256) :: dirname
629 dirname = trim(dir)//
'/hgh.'//trim(psp%atom_name)
630 call io_mkdir(trim(dirname), namespace)
631 hgh_unit =
io_open(trim(dirname)//
'/hgh', namespace, action=
'write')
632 loc_unit =
io_open(trim(dirname)//
'/local', namespace, action=
'write')
633 dat_unit =
io_open(trim(dirname)//
'/info', namespace, action=
'write')
634 kbp_unit =
io_open(trim(dirname)//
'/nonlocal', namespace, action=
'write')
635 wav_unit =
io_open(trim(dirname)//
'/wave', namespace, action=
'write')
638 write(hgh_unit,
'(a5,i6,5f12.6)') psp%atom_name, psp%z_val, psp%rlocal, psp%c(1:4)
639 write(hgh_unit,
'( 11x,4f12.6)') psp%rc(0), (psp%h(0,i,i), i = 1, 3)
641 write(hgh_unit,
'( 11x,4f12.6)') psp%rc(k), (psp%h(k, i, i), i = 1, 3)
642 write(hgh_unit,
'( 23x,4f12.6)') (psp%k(k, i, i), i = 1, 3)
646 write(dat_unit,
'(a,i3)')
'lmax = ', psp%l_max
647 if (psp%l_max >= 0)
then
648 write(dat_unit,
'(a,4f14.6)')
'kbr = ', psp%kbr
650 write(dat_unit,
'(a,5f14.6)')
'eigen = ', psp%eigen
651 write(dat_unit,
'(a,5f14.6)')
'occ = ', psp%conf%occ(1:psp%conf%p, 1)
653 do i = 1, psp%g%nrval
654 write(loc_unit, *) psp%g%rofi(i), psp%vlocal(i)
658 if (psp%l_max >= 0)
then
659 do i = 1, psp%g%nrval
660 write(kbp_unit,
'(10es14.4)') psp%g%rofi(i), ( (psp%kb(i, l, k) ,k = 1, 3),l = 0, psp%l_max)
665 do i = 1, psp%g%nrval
666 write(wav_unit, *) psp%g%rofi(i), (psp%rphi(i, l), l = 1, psp%conf%p)
Define which routines can be seen from the outside.
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public atomhxc(namespace, functl, g, nspin, dens, v, extra)
subroutine, public egofv(namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
egofv determines the eigenenergy and wavefunction corresponding ! to a particular l,...
real(real64), parameter, public m_two
real(real64), parameter, public m_max_exp_arg
real(real64), parameter, public m_huge
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_fourth
real(real64), parameter, public m_min_exp_arg
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
real(real64), parameter, public m_five
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public logrid_find_parameters(namespace, zz, aa, bb, np)
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 hgh_end(psp)
subroutine solve_schroedinger(psp, ierr, namespace)
subroutine projectorr_scalar(r, np, p, i, l, proj)
subroutine, public hgh_debug(psp, dir, namespace)
integer function load_params(unit, params, namespace)
subroutine, public hgh_get_eigen(psp, eigen)
subroutine, public hgh_process(psp, namespace)
subroutine, public hgh_init(psp, filename, namespace)
subroutine get_cutoff_radii(psp)
subroutine vlocalr_scalar(r, np, p, vloc)
The following data type contains: (a) the pseudopotential parameters, as read from a *....