26 use,
intrinsic :: iso_fortran_env
67 integer,
parameter,
public :: &
72 integer,
public,
parameter :: &
78 integer,
public,
parameter :: &
79 PROJ_J_INDEPENDENT = 0, & !< Non-relativisitic or scalar-relativistic pseudopotentials
83 integer,
parameter,
public :: INVALID_L = 333
85 character(len=4),
parameter :: ps_name(PSEUDO_FORMAT_UPF1:PSEUDO_FORMAT_PSP8) = &
86 (/
"upf1",
"upf2",
"qso ",
"psml",
"psf ",
"cpi ",
"fhi ",
"hgh ",
"psp8"/)
90 integer :: projector_type
91 integer :: relativistic_treatment
93 character(len=10),
private :: label
95 integer,
private :: ispin
96 real(real64),
private :: z
98 type(valconf_t) :: conf
99 type(logrid_t),
private :: g
100 type(spline_t),
allocatable :: ur(:, :)
101 type(spline_t),
allocatable,
private :: ur_sq(:, :)
102 logical,
allocatable :: bound(:, :)
109 logical :: no_vl = .false.
111 real(real64) :: projectors_sphere_threshold
118 real(real64) :: rc_max
121 integer :: projectors_per_l(5)
122 real(real64),
allocatable :: h(:,:,:), k(:, :, :)
123 type(spline_t),
allocatable :: kb(:, :)
124 type(spline_t),
allocatable :: dkb(:, :)
127 type(spline_t) :: core
128 type(spline_t) :: core_der
133 logical,
private :: has_long_range
135 type(spline_t),
private :: vlr
136 type(spline_t) :: vlr_sq
138 type(spline_t) :: nlr
140 real(real64) :: sigma_erf
142 logical,
private :: has_density
143 type(spline_t),
allocatable :: density(:)
144 type(spline_t),
allocatable :: density_der(:)
146 logical,
private :: is_separated
148 integer :: file_format
149 integer,
private :: pseudo_type
150 integer :: exchange_functional
151 integer :: correlation_functional
154 real(real64),
parameter :: eps = 1.0e-8_real64
160 subroutine ps_init(ps, namespace, label, z, user_lmax, user_llocal, ispin, filename)
161 type(
ps_t),
intent(out) :: ps
163 character(len=10),
intent(in) :: label
164 integer,
intent(in) :: user_lmax
165 integer,
intent(in) :: user_llocal
166 integer,
intent(in) :: ispin
167 real(real64),
intent(in) :: z
168 character(len=*),
intent(in) :: filename
170 integer :: l, ii, ll, is, ierr
176 real(real64),
allocatable :: eigen(:, :)
204 call parse_variable(namespace,
'SpeciesProjectorSphereThreshold', 0.001_real64, ps%projectors_sphere_threshold)
210 call messages_write(
"Cannot open pseudopotential file '"//trim(filename)//
"'.")
215 call messages_write(
"Cannot determine the pseudopotential type for species '"//trim(label)//
"' from", &
223 ps%relativistic_treatment = proj_j_independent
225 ps%sigma_erf = 0.625_real64
227 ps%projectors_per_l = 0
229 select case (ps%file_format)
233 ps%has_density = .false.
236 select case (ps%file_format)
247 ps%lmax = ps_psf%ps_grid%no_l_channels - 1
249 if (user_lmax /= invalid_l)
then
250 ps%lmax = min(ps%lmax, user_lmax)
251 if (user_lmax /= ps%lmax)
then
252 message(1) =
"lmax in Species block for " // trim(label) // &
253 " is larger than number available in pseudopotential."
258 ps%conf%p = ps_psf%ps_grid%no_l_channels
259 if (ps%lmax == 0) ps%llocal = 0
262 if (user_llocal == invalid_l)
then
265 ps%llocal = user_llocal
268 ps%projectors_per_l(1:ps%lmax+1) = 1
269 if (ps%llocal > -1) ps%projectors_per_l(ps%llocal+1) = 0
278 call ps_cpi_init(ps_cpi, trim(filename), namespace)
279 ps%conf%p = ps_cpi%ps_grid%no_l_channels
281 call ps_fhi_init(ps_fhi, trim(filename), namespace)
282 ps%conf%p = ps_fhi%ps_grid%no_l_channels
286 ps%conf%symbol = label(1:2)
295 ps%lmax = ps%conf%p - 1
297 if (user_lmax /= invalid_l)
then
298 ps%lmax = min(ps%lmax, user_lmax)
299 if (user_lmax /= ps%lmax)
then
300 message(1) =
"lmax in Species block for " // trim(label) // &
301 " is larger than number available in pseudopotential."
306 if (ps%lmax == 0) ps%llocal = 0
309 if (user_llocal == invalid_l)
then
312 ps%llocal = user_llocal
315 ps%projectors_per_l(1:ps%lmax+1) = 1
316 if (ps%llocal > -1) ps%projectors_per_l(ps%llocal+1) = 0
330 call hgh_init(ps_hgh, trim(filename), namespace)
334 ps%z_val = ps_hgh%z_val
337 ps%lmax = ps_hgh%l_max
338 ps%conf%symbol = label(1:2)
339 ps%sigma_erf = ps_hgh%rlocal
341 ps%projectors_per_l(1:ps%lmax+1) = 1
361 call ps_xml_init(ps_xml, namespace, trim(filename), ps%file_format, ierr)
370 if (ps_xml%kleinman_bylander)
then
371 ps%conf%p = ps_xml%nwavefunctions
373 ps%conf%p = ps_xml%lmax + 1
376 do ll = 0, ps_xml%lmax
377 ps%conf%l(ll + 1) = ll
380 ps%kbc = ps_xml%nchannels
381 ps%lmax = ps_xml%lmax
383 ps%projectors_per_l = 0
384 do ll = 0, ps_xml%lmax
388 if (ps_xml%kleinman_bylander)
then
389 ps%llocal = ps_xml%llocal
393 if (ps_xml%llocal >= 0) ps%llocal = ps_xml%llocal
394 if (user_llocal /= invalid_l) ps%llocal = user_llocal
395 assert(ps%llocal >= 0)
396 assert(ps%llocal <= ps%lmax)
399 ps%g%nrval = ps_xml%grid_size
401 safe_allocate(ps%g%rofi(1:ps%g%nrval))
402 safe_allocate(ps%g%r2ofi(1:ps%g%nrval))
404 do ii = 1, ps%g%nrval
405 ps%g%rofi(ii) = ps_xml%grid(ii)
406 ps%g%r2ofi(ii) = ps_xml%grid(ii)**2
411 ps%local = (ps%lmax == 0 .and. ps%llocal == 0 ) .or. (ps%lmax == -1 .and. ps%llocal == -1)
414 safe_allocate(ps%kb (0:ps%lmax, 1:ps%kbc))
415 safe_allocate(ps%dkb (0:ps%lmax, 1:ps%kbc))
416 safe_allocate(ps%ur (1:ps%conf%p, 1:ps%ispin))
417 safe_allocate(ps%ur_sq(1:ps%conf%p, 1:ps%ispin))
418 safe_allocate(ps%bound(1:ps%conf%p, 1:ps%ispin))
419 safe_allocate(ps%h (0:ps%lmax, 1:ps%kbc, 1:ps%kbc))
420 safe_allocate(ps%density(1:ps%ispin))
421 safe_allocate(ps%density_der(1:ps%ispin))
431 safe_allocate(eigen(1:ps%conf%p, 1:ps%ispin))
435 select case (ps%file_format)
448 safe_allocate(ps%k (0:ps%lmax, 1:ps%kbc, 1:ps%kbc))
458 call spline_der(ps%density(is), ps%density_der(is), ps%projectors_sphere_threshold)
463 call spline_der(ps%core, ps%core_der, ps%projectors_sphere_threshold)
468 ps%has_long_range = .
true.
469 ps%is_separated = .false.
473 safe_deallocate_a(eigen)
480 subroutine ps_info(ps, filename)
481 type(
ps_t),
intent(in) :: ps
482 character(len=*),
intent(in) :: filename
490 select case (ps%file_format)
491 case (pseudo_format_upf1)
499 case (pseudo_format_psp8)
521 select case (ps%pseudo_type)
535 select case (ps%file_format)
549 if (ps%llocal >= 0)
then
561 if (ps%llocal < 0)
then
590 type(
ps_t),
intent(inout) :: ps
592 real(real64),
allocatable :: vsr(:), vlr(:), nlr(:)
593 real(real64) :: r, exp_arg, arg, max_val_vsr
598 assert(ps%g%nrval > 0)
600 safe_allocate(vsr(1:ps%g%nrval))
601 safe_allocate(vlr(1:ps%g%nrval))
602 safe_allocate(nlr(1:ps%g%nrval))
608 do ii = 1, ps%g%nrval
611 if(arg > 5e-5_real64)
then
618 if(abs(vsr(ii)) < 1.0e-14_real64) vsr(ii) =
m_zero
619 max_val_vsr = max(max_val_vsr, abs(vsr(ii)))
621 exp_arg = -
m_half*r**2/ps%sigma_erf**2
630 call spline_fit(ps%g%nrval, ps%g%rofi, vlr, ps%vlr, ps%projectors_sphere_threshold)
633 call spline_fit(ps%g%nrval, ps%g%r2ofi, vlr, ps%vlr_sq, ps%projectors_sphere_threshold)
636 call spline_fit(ps%g%nrval, ps%g%rofi, nlr, ps%nlr, ps%projectors_sphere_threshold)
641 call spline_fit(ps%g%nrval, ps%g%rofi, vsr, ps%vl, ps%projectors_sphere_threshold)
642 if (max_val_vsr < 1.0e-12_real64) ps%no_vl = .
true.
644 safe_deallocate_a(vsr)
645 safe_deallocate_a(vlr)
646 safe_deallocate_a(nlr)
648 ps%is_separated = .
true.
656 type(
ps_t),
intent(inout) :: ps
664 if (l == ps%llocal) cycle
666 ps%rc_max = max(ps%rc_max, ps%kb(l, j)%x_threshold)
676 type(
ps_t),
intent(inout) :: ps
683 call spline_der(ps%kb(l, j), ps%dkb(l, j), ps%projectors_sphere_threshold)
693 type(
ps_t),
intent(inout) :: ps
694 integer,
intent(in) :: filter
695 real(real64),
intent(in) :: gmax
697 integer :: l, k, ispin
699 real(real64) :: alpha, beta_fs, rmax, rcut, gamma, beta_rs
705 case (ps_filter_none)
711 if(.not. ps%no_vl)
then
712 rmax = ps%vl%x_threshold
714 call spline_filter_mask(ps%vl, max(0, ps%llocal), rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
718 if (l == ps%llocal) cycle
720 call spline_filter_mask(ps%kb(l, k), l, ps%rc_max, gmax, alpha, gamma, ps%projectors_sphere_threshold)
725 rmax = ps%core%x_threshold
726 call spline_filter_mask(ps%core, 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
730 do ispin = 1, ps%ispin
732 rmax = ps%density(ispin)%x_threshold
733 call spline_filter_mask(ps%density(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
738 rmax = ps%density_der(ispin)%x_threshold
739 call spline_filter_mask(ps%density_der(ispin), 0, rmax, gmax, alpha, gamma, ps%projectors_sphere_threshold)
746 beta_fs = 18.0_real64
751 call spline_filter_bessel(ps%vl, max(0, ps%llocal), gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
753 if (l == ps%llocal) cycle
755 call spline_filter_bessel(ps%kb(l, k), l, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
760 call spline_filter_bessel(ps%core, 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
764 do ispin = 1, ps%ispin
765 call spline_filter_bessel(ps%density(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
767 call spline_filter_bessel(ps%density_der(ispin), 0, gmax, alpha, beta_fs, rcut, beta_rs, ps%projectors_sphere_threshold)
779 type(
ps_t),
intent(inout) :: ps
780 real(real64),
intent(in) :: eigen(:,:)
783 real(real64) :: ur1, ur2
797 if (.not. ps%bound(i, is)) cycle
799 do ir = ps%g%nrval, 3, -1
809 ur1 =
spline_eval(ps%ur(i, is), ps%g%rofi(ir-2))*ps%g%rofi(ir-2)
810 ur2 =
spline_eval(ps%ur(i, is), ps%g%rofi(ir-1))*ps%g%rofi(ir-1)
811 if ((ur1*ur2 >
m_zero) .and. (abs(ur2) > abs(ur1))) ps%bound(i, is) = .false.
823 subroutine ps_debug(ps, dir, namespace, gmax)
824 type(
ps_t),
intent(in) :: ps
825 character(len=*),
intent(in) :: dir
827 real(real64),
intent(in) :: gmax
830 type(
spline_t),
allocatable :: fw(:, :)
838 iunit =
io_open(trim(dir)//
'/pseudo-info', namespace, action=
'write')
839 write(iunit,
'(a,/)') ps%label
840 write(iunit,
'(a,a,/)')
'Format : ', ps_name(ps%file_format)
841 write(iunit,
'(a,f6.3)')
'z : ', ps%z
842 write(iunit,
'(a,f6.3,/)')
'zval : ', ps%z_val
843 write(iunit,
'(a,i4)')
'lmax : ', ps%lmax
844 write(iunit,
'(a,i4)')
'lloc : ', ps%llocal
845 write(iunit,
'(a,i4,/)')
'kbc : ', ps%kbc
846 write(iunit,
'(a,f9.5,/)')
'rcmax : ', ps%rc_max
847 write(iunit,
'(a,/)')
'h matrix:'
850 write(iunit,
'(10f9.5)') (ps%h(l, k, j), j = 1, ps%kbc)
853 if (
allocated(ps%k))
then
854 write(iunit,
'(/,a,/)')
'k matrix:'
857 write(iunit,
'(10f9.5)') (ps%k(l, k, j), j = 1, ps%kbc)
862 write(iunit,
'(/,a)')
'orbitals:'
864 if (ps%ispin == 2)
then
865 write(iunit,
'(1x,a,i2,3x,a,i2,3x,a,f5.1,3x,a,l1,3x,a,f5.1,3x,a,f5.1)')
'n = ', ps%conf%n(j),
'l = ', ps%conf%l(j), &
866 'j = ', ps%conf%j(j),
'bound = ', all(ps%bound(j,:)), &
867 'occ(1) = ', ps%conf%occ(j, 1),
'occ(2) = ', ps%conf%occ(j, 2)
869 write(iunit,
'(1x,a,i2,3x,a,i2,3x,a,f5.1,3x,a,l1,3x,a,f5.1)')
'n = ', ps%conf%n(j),
'l = ', ps%conf%l(j), &
870 'j = ', ps%conf%j(j),
'bound = ', all(ps%bound(j,:)), &
871 'occ = ', ps%conf%occ(j, 1)
879 iunit =
io_open(trim(dir)//
'/local', namespace, action=
'write')
884 iunit =
io_open(trim(dir)//
'/local_long_range', namespace, action=
'write')
889 iunit =
io_open(trim(dir)//
'/local_long_range_density', namespace, action=
'write')
894 iunit =
io_open(trim(dir)//
'/local_ft', namespace, action=
'write')
895 safe_allocate(fw(1, 1))
897 call spline_3dft(ps%vl, fw(1, 1), ps%projectors_sphere_threshold, gmax = gmax)
900 safe_deallocate_a(fw)
904 iunit =
io_open(trim(dir)//
'/nonlocal', namespace, action=
'write')
908 iunit =
io_open(trim(dir)//
'/nonlocal_derivative', namespace, action=
'write')
912 iunit =
io_open(trim(dir)//
'/nonlocal_ft', namespace, action=
'write')
913 safe_allocate(fw(0:ps%lmax, 1:ps%kbc))
917 call spline_besselft(ps%kb(k, j), fw(k, j), k, threshold=ps%projectors_sphere_threshold, &
923 safe_deallocate_a(fw)
927 iunit =
io_open(trim(dir)//
'/wavefunctions', namespace, action=
'write')
932 if (ps%has_density)
then
933 iunit =
io_open(trim(dir)//
'/density', namespace, action=
'write')
937 iunit =
io_open(trim(dir)//
'/density_derivative', namespace, action=
'write')
944 iunit =
io_open(trim(dir)//
'/nlcc', namespace, action=
'write')
955 type(
ps_t),
intent(inout) :: ps
957 if (.not.
allocated(ps%kb))
return
961 if (ps%is_separated)
then
976 if (
allocated(ps%density))
call spline_end(ps%density)
977 if (
allocated(ps%density_der))
call spline_end(ps%density_der)
981 safe_deallocate_a(ps%kb)
982 safe_deallocate_a(ps%dkb)
983 safe_deallocate_a(ps%ur)
984 safe_deallocate_a(ps%ur_sq)
985 safe_deallocate_a(ps%bound)
986 safe_deallocate_a(ps%h)
987 safe_deallocate_a(ps%k)
988 safe_deallocate_a(ps%density)
989 safe_deallocate_a(ps%density_der)
997 type(
ps_t),
intent(inout) :: ps
998 type(
ps_hgh_t),
intent(inout) :: ps_hgh
1004 if (ps%lmax >= 0)
then
1005 ps%rc_max = 1.1_real64 * maxval(ps_hgh%kbr(0:ps%lmax))
1009 ps%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%h(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1010 ps%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc) = ps_hgh%k(0:ps%lmax, 1:ps%kbc, 1:ps%kbc)
1021 integer :: l, is, j, ip
1022 real(real64),
allocatable :: hato(:), dens(:)
1026 safe_allocate(hato(1:ps_hgh%g%nrval))
1027 safe_allocate(dens(1:ps_hgh%g%nrval))
1030 do l = 0, ps_hgh%l_max
1032 hato = ps_hgh%kb(:, l, j)
1033 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%kb(l, j), ps%projectors_sphere_threshold)
1039 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, ps_hgh%vlocal, ps%vl, ps%projectors_sphere_threshold)
1046 hato(2:ps_hgh%g%nrval) = ps_hgh%rphi(2:ps_hgh%g%nrval, l)/ps_hgh%g%rofi(2:ps_hgh%g%nrval)
1049 do ip = 1, ps_hgh%g%nrval
1050 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(
m_four*
m_pi)
1053 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1054 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%r2ofi, hato, ps%ur_sq(l, is), ps%projectors_sphere_threshold)
1056 call spline_fit(ps_hgh%g%nrval, ps_hgh%g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1059 safe_deallocate_a(hato)
1060 safe_deallocate_a(dens)
1069 type(
ps_t),
intent(inout) :: ps
1075 ps%z_val = ps_grid%zval
1077 ps%nlcc = ps_grid%core_corrections
1079 ps%h(0:ps%lmax, 1, 1) = ps_grid%dkbcos(1:ps%lmax+1)
1083 ps%rc_max = maxval(ps_grid%kb_radius(1:ps%lmax+1)) * 1.5_real64
1089 ps%h(0:ps%lmax,:,:) = ps%h(0:ps%lmax,:,:) /
m_two
1098 real(real64),
allocatable :: hato(:), dens(:)
1099 integer :: is, l, ir, nrc, ip
1103 safe_allocate(hato(1:g%nrval))
1104 safe_allocate(dens(1:g%nrval))
1111 do l = 1, ps_grid%no_l_channels
1112 hato(2:) = ps_grid%rphi(2:, l, 1+is)/g%rofi(2:)
1117 dens(ip) = dens(ip) + ps%conf%occ(l, is)*hato(ip)**2/(
m_four*
m_pi)
1121 call spline_fit(g%nrval, g%rofi, hato, ps%ur(l, is), ps%projectors_sphere_threshold)
1122 call spline_fit(g%nrval, g%r2ofi, hato, ps%ur_sq(l, is), ps%projectors_sphere_threshold)
1125 call spline_fit(g%nrval, g%rofi, dens, ps%density(is), ps%projectors_sphere_threshold)
1132 hato(1:nrc) = ps_grid%KB(1:nrc, l)
1133 hato(nrc+1:g%nrval) =
m_zero
1135 call spline_fit(g%nrval, g%rofi, hato, ps%kb(l-1, 1), ps%projectors_sphere_threshold)
1140 hato(:) = ps_grid%vlocal(:)/
m_two
1141 call spline_fit(g%nrval, g%rofi, hato, ps%vl, ps%projectors_sphere_threshold)
1143 if (ps_grid%core_corrections)
then
1145 hato(2:) = ps_grid%chcore(2:)/(
m_four*
m_pi*g%r2ofi(2:))
1147 do ir = g%nrval-1, 2, -1
1148 if (hato(ir) > eps)
then
1154 hato(nrc:g%nrval) =
m_zero
1157 call spline_fit(g%nrval, g%rofi, hato, ps%core, ps%projectors_sphere_threshold)
1160 safe_deallocate_a(hato)
1161 safe_deallocate_a(dens)
1170 type(
ps_t),
intent(inout) :: ps
1171 type(
ps_xml_t),
intent(in) :: ps_xml
1174 integer :: ll, ip, is, ic, jc, ir, nrc, ii
1175 real(real64) :: rr, kbcos, kbnorm, dnrm, avgv, volume_element
1176 real(real64),
allocatable :: vlocal(:), kbprojector(:), wavefunction(:), nlcc_density(:), dens(:,:)
1177 integer,
allocatable :: cmap(:, :)
1178 real(real64),
allocatable :: matrix(:, :), eigenvalues(:)
1179 logical :: density_is_known
1180 integer,
allocatable :: order(:)
1181 real(real64),
allocatable :: occ_tmp(:,:)
1188 if (ps%file_format == pseudo_format_psp8)
then
1190 message(1) =
"SOC from PSP8 is not currently supported."
1191 message(2) =
"Only scalar relativistic effects will be considered."
1198 ps%nlcc = ps_xml%nlcc
1200 ps%z_val = ps_xml%valence_charge
1202 ps%has_density = ps_xml%has_density
1205 safe_allocate(vlocal(1:ps%g%nrval))
1207 do ip = 1, ps%g%nrval
1208 rr = ps_xml%grid(ip)
1209 if (ip <= ps_xml%grid_size)
then
1210 vlocal(ip) = ps_xml%potential(ip, ps%llocal)
1212 vlocal(ip) = -ps_xml%valence_charge/rr
1216 call spline_fit(ps%g%nrval, ps%g%rofi, vlocal, ps%vl, ps%projectors_sphere_threshold)
1218 safe_deallocate_a(vlocal)
1220 safe_allocate(kbprojector(1:ps%g%nrval))
1221 safe_allocate(wavefunction(1:ps%g%nrval))
1226 density_is_known = .false.
1229 if (ps_xml%kleinman_bylander)
then
1231 safe_allocate(cmap(0:ps_xml%lmax, 1:ps_xml%nchannels))
1234 do ll = 0, ps_xml%lmax
1235 do ic = 1, ps_xml%nchannels
1239 if (ll == ps_xml%llocal) cycle
1243 assert(mod(ps_xml%nchannels, 2)==0)
1245 cmap(ll, ic) = int((ic-1)/2)*2 + 2
1248 cmap(ll, ic) = int((ic-1)/2)*2 + 1
1254 assert(sum(cmap(ll, 1:ps_xml%nchannels)) == (ps_xml%nchannels + 1)*ps_xml%nchannels/2)
1257 assert(all(cmap >= 0 .and. cmap <= ps_xml%nchannels))
1259 safe_allocate(matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels))
1260 safe_allocate(eigenvalues(1:ps_xml%nchannels))
1265 do ll = 0, ps_xml%lmax
1267 if (
is_diagonal(ps_xml%nchannels, ps_xml%dij(ll, :, :)) .or. &
1270 do ic = 1, ps_xml%nchannels
1271 eigenvalues(ic) = ps_xml%dij(ll, ic, ic)
1272 matrix(ic, ic) =
m_one
1276 matrix(1:ps_xml%nchannels, 1:ps_xml%nchannels) = ps_xml%dij(ll, 1:ps_xml%nchannels, 1:ps_xml%nchannels)
1280 do ic = 1, ps_xml%nchannels
1282 do ip = 1, ps%g%nrval
1284 if (ip <= ps_xml%grid_size)
then
1285 do jc = 1, ps_xml%nchannels
1286 kbprojector(ip) = kbprojector(ip) + matrix(jc, ic)*ps_xml%projector(ip, ll, jc)
1291 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, cmap(ll, ic)), ps%projectors_sphere_threshold)
1293 ps%h(ll, cmap(ll, ic), cmap(ll, ic)) = eigenvalues(ic)
1299 safe_deallocate_a(matrix)
1300 safe_deallocate_a(eigenvalues)
1302 ps%conf%p = ps_xml%nwavefunctions
1309 if ((.not.
ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0)
then
1310 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1314 do ii = 1, ps_xml%nwavefunctions
1316 ps%conf%n(ii) = ps_xml%wf_n(ii)
1317 ps%conf%l(ii) = ps_xml%wf_l(ii)
1319 if (ps%ispin == 2)
then
1320 ps%conf%occ(ii, 1) = min(ps_xml%wf_occ(ii),
m_two*ps_xml%wf_l(ii) +
m_one)
1321 ps%conf%occ(ii, 2) = ps_xml%wf_occ(ii) - ps%conf%occ(ii, 1)
1323 ps%conf%occ(ii, 1) = ps_xml%wf_occ(ii)
1331 do ip = 1, ps%g%nrval
1332 if (ip <= ps_xml%grid_size)
then
1333 wavefunction(ip) = ps_xml%wavefunction(ip, ii)
1335 wavefunction(ip) =
m_zero
1340 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ii, is), ps%projectors_sphere_threshold)
1341 call spline_fit(ps%g%nrval, ps%g%r2ofi, wavefunction, ps%ur_sq(ii, is), ps%projectors_sphere_threshold)
1346 do ip = 1, ps_xml%grid_size
1347 dens(ip, is) = dens(ip, is) + ps%conf%occ(ii, is)*wavefunction(ip)**2/(
m_four*
m_pi)
1355 if ((.not.
ps_has_density(ps) .or. ps%ispin == 2) .and. ps_xml%nwavefunctions > 0)
then
1357 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1359 safe_deallocate_a(dens)
1360 density_is_known = .
true.
1361 ps%has_density = .
true.
1364 safe_deallocate_a(cmap)
1370 ps%conf%symbol = ps%label(1:2)
1374 safe_allocate(order(1:ps_xml%lmax+1))
1375 safe_allocate(occ_tmp(1:ps_xml%lmax+1, 1:2))
1376 occ_tmp = ps%conf%occ
1377 call sort(ps%conf%l(1:ps_xml%lmax+1), order)
1378 do ll = 0, ps_xml%lmax
1379 ps%conf%occ(ll+1, 1:2) = occ_tmp(order(ll+1), 1:2)
1381 safe_deallocate_a(order)
1382 safe_deallocate_a(occ_tmp)
1385 if (abs(sum(ps%conf%occ) - ps%z_val ) <
m_epsilon)
then
1386 safe_allocate(dens(1:ps%g%nrval, 1:ps%ispin))
1390 do ll = 0, ps_xml%lmax
1395 do ip = 1, ps_xml%grid_size
1396 rr = ps_xml%grid(ip)
1397 volume_element = rr**2*ps_xml%weights(ip)
1398 kbprojector(ip) = (ps_xml%potential(ip, ll) - ps_xml%potential(ip, ps%llocal))*ps_xml%wavefunction(ip, ll)
1399 dnrm = dnrm + kbprojector(ip)**2*volume_element
1400 avgv = avgv + kbprojector(ip)*ps_xml%wavefunction(ip, ll)*volume_element
1403 kbcos = dnrm/(safe_tol(avgv,1.0e-20_real64))
1404 kbnorm =
m_one/(safe_tol(
sqrt(dnrm),1.0e-20_real64))
1406 if (ll /= ps%llocal)
then
1407 ps%h(ll, 1, 1) = kbcos
1408 kbprojector = kbprojector*kbnorm
1413 call spline_fit(ps%g%nrval, ps%g%rofi, kbprojector, ps%kb(ll, 1), ps%projectors_sphere_threshold)
1416 do ip = 1, ps%g%nrval
1417 if (ip <= ps_xml%grid_size)
then
1418 wavefunction(ip) = ps_xml%wavefunction(ip, ll)
1420 wavefunction(ip) =
m_zero
1425 call spline_fit(ps%g%nrval, ps%g%rofi, wavefunction, ps%ur(ll + 1, is), ps%projectors_sphere_threshold)
1426 call spline_fit(ps%g%nrval, ps%g%r2ofi, wavefunction, ps%ur_sq(ll + 1, is), ps%projectors_sphere_threshold)
1430 if (abs(sum(ps%conf%occ) - ps%z_val) <
m_epsilon)
then
1432 do ip = 1, ps_xml%grid_size
1433 dens(ip, is) = dens(ip, is) + ps%conf%occ(ll+1, is)*wavefunction(ip)**2/(
m_four*
m_pi)
1440 if (abs(sum(ps%conf%occ) - ps%z_val) <
m_epsilon)
then
1442 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:,is), ps%density(is), ps%projectors_sphere_threshold)
1444 safe_deallocate_a(dens)
1445 density_is_known = .
true.
1446 ps%has_density = .
true.
1453 safe_allocate(dens(1:ps%g%nrval, 1))
1455 dens(1:ps_xml%grid_size, 1) = ps_xml%density(1:ps_xml%grid_size)/ps%ispin
1456 dens(ps_xml%grid_size + 1:ps%g%nrval, 1) =
m_zero
1459 call spline_fit(ps%g%nrval, ps%g%rofi, dens(:, 1), ps%density(is), ps%projectors_sphere_threshold)
1462 safe_deallocate_a(dens)
1466 if (ps_xml%nlcc)
then
1468 safe_allocate(nlcc_density(1:ps%g%nrval))
1470 nlcc_density(1:ps_xml%grid_size) = ps_xml%nlcc_density(1:ps_xml%grid_size)
1473 do ir = ps_xml%grid_size - 1, 1, -1
1474 if (nlcc_density(ir) > eps)
then
1480 nlcc_density(nrc:ps%g%nrval) =
m_zero
1482 call spline_fit(ps%g%nrval, ps%g%rofi, nlcc_density, ps%core, ps%projectors_sphere_threshold)
1484 safe_deallocate_a(nlcc_density)
1490 ps%rc_max = ps%rc_max*1.05_real64
1492 safe_deallocate_a(kbprojector)
1493 safe_deallocate_a(wavefunction)
1501 integer,
intent(in) :: dim
1502 real(real64),
intent(in) :: matrix(:, :)
1510 if (abs(matrix(ii, jj)) > 1e10_real64)
is_diagonal = .false.
1519 type(
ps_t),
intent(in) :: ps
1534 type(
ps_t),
intent(in) :: ps
1541 if (any(.not. ps%bound(i,:))) cycle
1550 type(
ps_t),
intent(in) :: ps
1552 has_density = ps%has_density
1558 pure logical function ps_has_nlcc(ps)
result(has_nlcc)
1559 type(
ps_t),
intent(in) :: ps
1566 real(real64) function ps_density_volume(ps, namespace) result(volume)
1567 type(
ps_t),
intent(in) :: ps
1570 integer :: ip, ispin
1572 real(real64),
allocatable ::vol(:)
1575 push_sub(ps_density_volume)
1578 message(1) =
"The pseudopotential does not contain an atomic density"
1582 safe_allocate(vol(1:ps%g%nrval))
1584 do ip = 1, ps%g%nrval
1587 do ispin = 1, ps%ispin
1593 call spline_fit(ps%g%nrval, ps%g%rofi, vol, volspl, ps%projectors_sphere_threshold)
1597 safe_deallocate_a(vol)
1599 pop_sub(ps_density_volume)
1606 type(namespace_t),
intent(in) :: namespace
1607 real(real64),
intent(in) :: zz
1608 real(real64),
intent(in) :: valcharge
1609 integer,
intent(in) :: ispin
1610 type(valconf_t),
intent(inout) :: conf
1622 if (debug%info)
then
1623 write(message(1),
'(a,a,a)')
'Debug: Guessing the atomic occupations for ', trim(conf%symbol),
"."
1624 call messages_info(1, namespace=namespace)
1627 assert(valcharge <= zz)
1631 if(int(zz) > 2 .and. val > zz - 2)
then
1635 if(int(zz) > 4 .and. val > zz - 4)
then
1640 if(int(zz) > 18 .and. val > zz - 10)
then
1644 if(int(zz) > 12 .and. val > zz - 12)
then
1648 if(int(zz) > 18 .and. val > zz - 18)
then
1652 if(int(zz) > 28 .and. val > zz - 28)
then
1656 if(int(zz) > 30 .and. val > zz - 30)
then
1660 if(int(zz) > 36 .and. val > zz - 36)
then
1664 if(int(zz) > 46 .and. val > zz - 46)
then
1669 if((int(zz) > 48 .and. val > zz - 48) .or. &
1670 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 62))
then
1675 if((int(zz) > 54 .and. val > zz - 54) .or. &
1676 (int(zz) > 70 .and. int(zz) <= 81 .and. val > zz - 68) )
then
1681 if(int(zz) > 80 .and. val > zz - 68)
then
1687 select case (int(zz))
1934 if (val < m_epsilon)
then
1936 if (ispin == 2)
then
1937 call valconf_unpolarized_to_polarized(conf)
1941 message(1) =
"Error in attributing atomic occupations"
1942 call messages_warning(1, namespace=namespace)
1949 real(real64),
intent(inout) :: val
1950 integer,
intent(in) :: max_occ, nn
1953 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
1954 val = val - conf%occ(conf%p,1)
1960 real(real64),
intent(inout) :: val
1961 integer,
intent(in) :: max_occ, nn
1964 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
1965 val = val - conf%occ(conf%p,1)
1971 real(real64),
intent(inout) :: val
1972 integer,
intent(in) :: max_occ, nn
1975 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
1976 val = val - conf%occ(conf%p,1)
1982 real(real64),
intent(inout) :: val
1983 integer,
intent(in) :: max_occ, nn
1986 conf%occ(conf%p,1) = min(val, real(max_occ, real64) )
1987 val = val - conf%occ(conf%p,1)
1994#include "ps_pspio_inc.F90"
This is the common interface to a sorting routine. It performs the shell algorithm,...
Some operations may be done for one spline-function, or for an array of them.
The integral may be done with or without integration limits, but we want the interface to be common.
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public valconf_unpolarized_to_polarized(conf)
subroutine, public valconf_copy(cout, cin)
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_min_exp_arg
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
subroutine, public io_close(iunit, grp)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public logrid_copy(grid_in, grid_out)
integer function, public logrid_index(grid, rofi)
subroutine, public logrid_end(grid)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
subroutine, public messages_new_line()
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 messages_input_error(namespace, var, details, row, column)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public ps_cpi_end(ps_cpi)
subroutine, public ps_cpi_init(ps_cpi, filename, namespace)
subroutine, public ps_cpi_process(ps_cpi, lloc, namespace)
subroutine, public ps_fhi_end(ps_fhi)
subroutine, public ps_fhi_init(ps_fhi, filename, namespace)
subroutine, public ps_fhi_process(ps_fhi, lmax, lloc, namespace)
subroutine, public hgh_end(psp)
subroutine, public hgh_get_eigen(psp, eigen)
subroutine, public hgh_process(psp, namespace)
subroutine, public hgh_init(psp, filename, namespace)
real(real64) function, public first_point_extrapolate(x, y, high_order)
pure logical function, public ps_has_density(ps)
integer, parameter, public ps_filter_ts
subroutine, public ps_end(ps)
subroutine hgh_load(ps, ps_hgh)
subroutine, public ps_init(ps, namespace, label, z, user_lmax, user_llocal, ispin, filename)
subroutine, public ps_separate(ps)
separate the local potential into (soft) long-ranged and (hard) short-ranged parts
subroutine ps_check_bound(ps, eigen)
subroutine ps_grid_load(ps, ps_grid)
subroutine, public ps_debug(ps, dir, namespace, gmax)
integer, parameter, public proj_hgh
pure integer function, public ps_bound_niwfs(ps)
Returns the number of bound atomic orbitals taking into account then m quantum number multiplicity.
subroutine, public ps_getradius(ps)
subroutine, public ps_guess_atomic_occupations(namespace, zz, valcharge, ispin, conf)
This routines provides, given Z and the number of valence electron the occupations of the orbitals....
pure logical function, public ps_has_nlcc(ps)
pure integer function, public ps_niwfs(ps)
Returns the number of atomic orbitals taking into account then m quantum number multiplicity.
subroutine, public ps_pspio_init(ps, namespace, label, z, lmax, ispin, filename)
subroutine ps_info(ps, filename)
integer, parameter, public proj_j_average
Fully-relativistic pseudopotentials with separate j-average and SOC terms.
integer, parameter, public proj_rkb
integer, parameter, public proj_j_dependent
Fully-relativistic j-dependent pseudopotentials.
subroutine, public ps_derivatives(ps)
subroutine, public ps_filter(ps, filter, gmax)
real(real64) function, public ps_density_volume(ps, namespace)
integer, parameter, public ps_filter_bsb
integer, parameter, public proj_kb
subroutine ps_xml_load(ps, ps_xml, namespace)
logical function is_diagonal(dim, matrix)
subroutine, public ps_psf_init(pstm, ispin, filename, namespace)
subroutine, public ps_psf_process(ps_psf, namespace, lmax, lloc)
subroutine, public ps_psf_get_eigen(ps_psf, eigen)
subroutine, public ps_psf_end(ps_psf)
subroutine, public ps_xml_end(this)
subroutine, public ps_xml_init(this, namespace, filename, fmt, ierr)
integer, parameter, public pseudo_format_file_not_found
integer, parameter, public pseudo_type_kleinman_bylander
integer, parameter, public pseudo_format_qso
integer, parameter, public pseudo_format_hgh
integer, parameter, public pseudo_format_upf2
integer, parameter, public pseudo_type_paw
integer, parameter, public pseudo_format_cpi
integer, parameter, public pseudo_exchange_unknown
integer, parameter, public pseudo_type_semilocal
integer, parameter, public pseudo_type_ultrasoft
logical pure function, public pseudo_has_total_angular_momentum(pseudo)
integer, parameter, public pseudo_correlation_unknown
integer, parameter, public pseudo_format_psf
integer, parameter, public pseudo_format_fhi
integer, parameter, public pseudo_format_unknown
integer, parameter, public pseudo_format_psml
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public spline_filter_mask(spl, l, rmax, qmax, alpha, gamma, threshold)
subroutine, public spline_filter_bessel(spl, l, qmax, alpha, beta_fs, rcut, beta_rs, threshold)
subroutine, public spline_der(spl, dspl, threshold)
subroutine, public spline_force_pos(spl, threshold)
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
subroutine, public spline_3dft(spl, splw, threshold, gmax)
subroutine, public spline_besselft(spl, splw, l, threshold, gmax)
real(real64) function, public spline_eval(spl, x)
subroutine fill_s_orbs(val, max_occ, nn)
subroutine fill_p_orbs(val, max_occ, nn)
subroutine fill_d_orbs(val, max_occ, nn)
subroutine fill_f_orbs(val, max_occ, nn)
remember that the FHI format is basically the CPI format with a header
The following data type contains: (a) the pseudopotential parameters, as read from a *....
the basic spline datatype