38 use,
intrinsic :: iso_fortran_env
101 logical :: calculating
102 logical :: time_present
104 real(real64),
allocatable :: density(:, :)
105 logical :: total_density_alloc
106 real(real64),
pointer,
contiguous :: total_density(:)
107 type(energy_t),
allocatable :: energy
108 type(states_elec_t),
pointer :: hf_st
109 real(real64),
allocatable :: vxc(:, :)
110 real(real64),
allocatable :: vtau(:, :)
111 real(real64),
allocatable :: axc(:, :, :)
112 real(real64),
allocatable :: a_ind(:, :)
113 real(real64),
allocatable :: b_ind(:, :)
114 logical :: calc_energy
119 integer,
public :: theory_level = -1
121 logical,
public :: frozen_hxc = .false.
123 integer,
public :: xc_family = 0
124 integer,
public :: xc_flags = 0
125 integer,
public :: xc_photon = 0
126 type(xc_t),
public :: xc
127 type(xc_photons_t),
public :: xc_photons
128 type(xc_oep_t),
public :: oep
129 type(xc_oep_photon_t),
public :: oep_photon
130 type(xc_ks_inversion_t),
public :: ks_inversion
131 type(xc_sic_t),
public :: sic
132 type(xc_vdw_t),
public :: vdw
133 type(grid_t),
pointer,
public :: gr
134 type(v_ks_calc_t) :: calc
135 logical :: calculate_current = .false.
136 type(current_t) :: current_calculator
137 logical :: include_td_field = .false.
138 logical,
public :: has_photons = .false.
139 logical :: xc_photon_include_hartree = .
true.
141 real(real64),
public :: stress_xc_gga(3, 3)
142 type(photon_mode_t),
pointer,
public :: pt => null()
143 type(mf_t),
public :: pt_mx
149 subroutine v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
150 type(v_ks_t),
intent(inout) :: ks
151 type(namespace_t),
intent(in) :: namespace
152 type(grid_t),
target,
intent(inout) :: gr
153 type(states_elec_t),
intent(in) :: st
154 type(ions_t),
intent(inout) :: ions
155 type(multicomm_t),
intent(in) :: mc
156 class(space_t),
intent(in) :: space
157 type(kpoints_t),
intent(in) :: kpoints
159 integer :: x_id, c_id, xk_id, ck_id, default, val
160 logical :: parsed_theory_level, using_hartree_fock
161 integer :: pseudo_x_functional, pseudo_c_functional
204 ks%xc_family = xc_family_none
210 parsed_theory_level = .false.
231 call messages_write(
'Info: the XCFunctional has been selected to match the pseudopotentials', new_line = .
true.)
246 call messages_write(
'The XCFunctional that you selected does not match the one used', new_line = .
true.)
295 call parse_variable(namespace,
'XCPhotonFunctional', option__xcphotonfunctional__none, ks%xc_photon)
305 call parse_variable(namespace,
'XCPhotonIncludeHartree', .
true., ks%xc_photon_include_hartree)
307 if (.not. ks%xc_photon_include_hartree)
then
318 using_hartree_fock = (ks%theory_level ==
hartree_fock) &
320 call xc_init(ks%xc, namespace, space%dim, space%periodic_dim, st%qtot, &
321 x_id, c_id, xk_id, ck_id,
hartree_fock = using_hartree_fock, ispin=st%d%ispin)
323 ks%xc_family = ks%xc%family
324 ks%xc_flags = ks%xc%flags
326 if (.not. parsed_theory_level)
then
335 call parse_variable(namespace,
'TheoryLevel', default, ks%theory_level)
345 if (
accel_is_enabled() .and. (gr%parallel_in_domains .or. st%parallel_in_states .or. st%d%kpt%parallel))
then
352 ks%xc_family = ior(ks%xc_family, xc_family_oep)
362 ks%sic%amaldi_factor =
m_one
364 select case (ks%theory_level)
369 if (space%periodic_dim == space%dim)
then
372 if (kpoints%full%npoints > 1)
then
377 if (kpoints%full%npoints > 1)
then
392 if (
bitand(ks%xc_family, xc_family_lda + xc_family_gga) /= 0)
then
393 call xc_sic_init(ks%sic, namespace, gr, st, mc, space)
396 if (
bitand(ks%xc_family, xc_family_oep) /= 0)
then
397 select case (ks%xc%functional(
func_x,1)%id)
399 if (kpoints%reduced%npoints > 1)
then
404 if (kpoints%reduced%npoints > 1)
then
409 if((.not. ks%has_photons) .or. (ks%xc_photon /= 0))
then
410 if(oep_type == -1)
then
413 call xc_oep_init(ks%oep, namespace, gr, st, mc, space, oep_type)
427 message(1) =
"SICCorrection can only be used with Kohn-Sham DFT"
431 if (st%d%ispin ==
spinors)
then
432 if (
bitand(ks%xc_family, xc_family_mgga + xc_family_hyb_mgga) /= 0)
then
437 ks%frozen_hxc = .false.
442 ks%calc%calculating = .false.
447 call ks%vdw%init(namespace, space, gr, ks%xc, ions, x_id, c_id)
449 if (ks%xc_photon /= 0)
then
451 call ks%xc_photons%init(namespace, ks%xc_photon , space, gr, st)
462 integer,
intent(out) :: x_functional
463 integer,
intent(out) :: c_functional
465 integer :: xf, cf, ispecies
466 logical :: warned_inconsistent
471 warned_inconsistent = .false.
472 do ispecies = 1, ions%nspecies
473 select type(spec=>ions%species(ispecies)%s)
475 xf = spec%x_functional()
476 cf = spec%c_functional()
483 call messages_write(
"Unknown XC functional for species '"//trim(ions%species(ispecies)%s%get_label())//
"'")
491 if (xf /= x_functional .and. .not. warned_inconsistent)
then
492 call messages_write(
'Inconsistent XC functional detected between species')
494 warned_inconsistent = .
true.
501 if (cf /= c_functional .and. .not. warned_inconsistent)
then
502 call messages_write(
'Inconsistent XC functional detected between species')
504 warned_inconsistent = .
true.
519 type(
v_ks_t),
intent(inout) :: ks
525 select case (ks%theory_level)
530 if (
bitand(ks%xc_family, xc_family_oep) /= 0)
then
540 if (ks%xc_photon /= 0)
then
541 call ks%xc_photons%end()
551 type(
v_ks_t),
intent(in) :: ks
552 integer,
optional,
intent(in) :: iunit
553 type(
namespace_t),
optional,
intent(in) :: namespace
560 select case (ks%theory_level)
585 subroutine v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
588 type(
grid_t),
intent(in) :: gr
589 type(
ions_t),
intent(in) :: ions
592 type(
v_ks_t),
intent(inout) :: ks
594 logical,
optional,
intent(in) :: calc_eigenval
595 logical,
optional,
intent(in) :: calc_current
597 integer,
allocatable :: ind(:)
599 real(real64),
allocatable :: copy_occ(:)
600 logical :: calc_eigenval_
601 logical :: calc_current_
609 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = calc_eigenval_, calc_current = calc_current_)
611 if (st%restart_reorder_occs .and. .not. st%fromScratch)
then
612 message(1) =
"Reordering occupations for restart."
615 safe_allocate(ind(1:st%nst))
616 safe_allocate(copy_occ(1:st%nst))
619 call sort(st%eigenval(:, ik), ind)
620 copy_occ(1:st%nst) = st%occ(1:st%nst, ik)
622 st%occ(ist, ik) = copy_occ(ind(ist))
626 safe_deallocate_a(ind)
627 safe_deallocate_a(copy_occ)
637 subroutine v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
638 calc_eigenval, time, calc_energy, calc_current, force_semilocal)
639 type(
v_ks_t),
intent(inout) :: ks
644 type(
ions_t),
intent(in) :: ions
646 logical,
optional,
intent(in) :: calc_eigenval
647 real(real64),
optional,
intent(in) :: time
648 logical,
optional,
intent(in) :: calc_energy
649 logical,
optional,
intent(in) :: calc_current
650 logical,
optional,
intent(in) :: force_semilocal
652 logical :: calc_current_
658 call v_ks_calc_start(ks, namespace, space, hm, st, ions, hm%kpoints%latt, ext_partners, time, &
659 calc_energy, calc_current_, force_semilocal=force_semilocal)
661 ext_partners, force_semilocal=force_semilocal)
679 subroutine v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, &
680 calc_energy, calc_current, force_semilocal)
681 type(
v_ks_t),
target,
intent(inout) :: ks
683 class(
space_t),
intent(in) :: space
686 type(
ions_t),
intent(in) :: ions
689 real(real64),
optional,
intent(in) :: time
690 logical,
optional,
intent(in) :: calc_energy
691 logical,
optional,
intent(in) :: calc_current
692 logical,
optional,
intent(in) :: force_semilocal
694 logical :: calc_current_
699 .and. ks%calculate_current &
705 assert(.not. ks%calc%calculating)
706 ks%calc%calculating = .
true.
709 write(
message(1),
'(a)')
'Debug: Calculating Kohn-Sham potential.'
713 ks%calc%time_present =
present(time)
719 if (ks%frozen_hxc)
then
720 if (calc_current_)
then
730 safe_allocate(ks%calc%energy)
734 ks%calc%energy%intnvxc =
m_zero
736 nullify(ks%calc%total_density)
746 if (ks%theory_level /=
hartree .and. ks%theory_level /=
rdmft)
call v_a_xc(hm, force_semilocal)
748 ks%calc%total_density_alloc = .false.
751 if (calc_current_)
then
756 nullify(ks%calc%hf_st)
760 safe_allocate(ks%calc%hf_st)
763 if (st%parallel_in_states)
then
765 call messages_write(
'State parallelization of Hartree-Fock exchange is not supported')
767 call messages_write(
'when running with OpenCL/CUDA. Please use domain parallelization')
769 call messages_write(
"or disable acceleration using 'DisableAccel = yes'.")
781 if (hm%self_induced_magnetic)
then
782 safe_allocate(ks%calc%a_ind(1:ks%gr%np_part, 1:space%dim))
783 safe_allocate(ks%calc%b_ind(1:ks%gr%np_part, 1:space%dim))
784 call magnetic_induced(namespace, ks%gr, st, hm%psolver, hm%kpoints, ks%calc%a_ind, ks%calc%b_ind)
787 if ((ks%has_photons) .and. (ks%calc%time_present) .and. (ks%xc_photon == 0) )
then
788 call mf_calc(ks%pt_mx, ks%gr, st, ions, ks%pt, time)
806 safe_allocate(ks%calc%density(1:ks%gr%np, 1:st%d%nspin))
811 call lalg_scal(ks%gr%np, st%d%nspin, ks%sic%amaldi_factor, ks%calc%density)
814 nullify(ks%calc%total_density)
815 if (
allocated(st%rho_core) .or. hm%d%spin_channels > 1)
then
816 ks%calc%total_density_alloc = .
true.
818 safe_allocate(ks%calc%total_density(1:ks%gr%np))
821 ks%calc%total_density(ip) = sum(ks%calc%density(ip, 1:hm%d%spin_channels))
825 if (
allocated(st%rho_core))
then
826 call lalg_axpy(ks%gr%np, -ks%sic%amaldi_factor, st%rho_core, ks%calc%total_density)
829 ks%calc%total_density_alloc = .false.
830 ks%calc%total_density => ks%calc%density(:, 1)
837 subroutine v_a_xc(hm, force_semilocal)
839 logical,
optional,
intent(in) :: force_semilocal
846 ks%calc%energy%exchange =
m_zero
847 ks%calc%energy%correlation =
m_zero
848 ks%calc%energy%xc_j =
m_zero
849 ks%calc%energy%vdw =
m_zero
851 safe_allocate(ks%calc%vxc(1:ks%gr%np, 1:st%d%nspin))
855 safe_allocate(ks%calc%vtau(1:ks%gr%np, 1:st%d%nspin))
860 if (ks%calc%calc_energy)
then
862 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, st%d%ispin, &
863 latt%rcell_volume, ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation, &
864 deltaxc = ks%calc%energy%delta_xc, vtau = ks%calc%vtau, force_orbitalfree=force_semilocal)
866 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, st%d%ispin, &
867 latt%rcell_volume, ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation, &
868 deltaxc = ks%calc%energy%delta_xc, stress_xc=ks%stress_xc_gga, force_orbitalfree=force_semilocal)
872 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, &
873 st%d%ispin, latt%rcell_volume, ks%calc%vxc, vtau = ks%calc%vtau, force_orbitalfree=force_semilocal)
875 call xc_get_vxc(ks%gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, ks%calc%density, &
876 st%d%ispin, latt%rcell_volume, ks%calc%vxc, stress_xc=ks%stress_xc_gga, force_orbitalfree=force_semilocal)
882 if (st%d%ispin /=
spinors)
then
883 message(1) =
"Noncollinear functionals can only be used with spinor wavefunctions."
888 message(1) =
"Cannot perform LCAO for noncollinear MGGAs."
889 message(2) =
"Please perform a LDA calculation first."
893 if (ks%calc%calc_energy)
then
895 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc, &
896 vtau = ks%calc%vtau, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
898 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc, &
899 ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
903 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, &
904 ks%calc%vxc, vtau = ks%calc%vtau)
906 call xc_get_nc_vxc(ks%gr, ks%xc, st, hm%kpoints, space, namespace, ks%calc%density, ks%calc%vxc)
922 if (ks%calc%calc_energy)
then
923 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
924 ks%calc%vxc, ex = ks%calc%energy%exchange, ec = ks%calc%energy%correlation)
926 call xc_sic_calc_adsic(ks%sic, namespace, space, ks%gr, st, hm, ks%xc, ks%calc%density, &
938 call x_slater_calc(namespace, ks%gr, space, hm%exxop, st, hm%kpoints, ks%calc%energy%exchange, &
941 call x_fbe_calc(ks%xc%functional(
func_x,1)%id, namespace, hm%psolver, ks%gr, st, space, &
942 ks%calc%energy%exchange, vxc = ks%calc%vxc)
946 call fbe_c_lda_sl(namespace, hm%psolver, ks%gr, st, space, &
947 ks%calc%energy%correlation, vxc = ks%calc%vxc)
953 hm, st, space, ks%calc%energy%exchange, ks%calc%energy%correlation, vxc = ks%calc%vxc)
956 hm, st, space, ks%calc%energy%exchange, ks%calc%energy%correlation, vxc = ks%calc%vxc)
958 ks%calc%energy%photon_exchange = ks%oep_photon%pt%ex
966 call xc_ks_inversion_calc(ks%ks_inversion, namespace, space, ks%gr, hm, ext_partners, st, vxc = ks%calc%vxc, &
971 if (ks%xc_photon /= 0)
then
973 call ks%xc_photons%v_ks(namespace, ks%calc%total_density, ks%calc%density, ks%gr, space, hm%psolver, hm%ep, st)
976 do ispin = 1, hm%d%spin_channels
977 call lalg_axpy(ks%gr%np,
m_one, ks%xc_photons%vpx(1:ks%gr%np), ks%calc%vxc(1:ks%gr%np, ispin) )
981 ks%calc%energy%photon_exchange = ks%xc_photons%ex
986 call ks%vdw%calc(namespace, space, latt, ions%atom, ions%natoms, ions%pos, &
987 ks%gr, st, ks%calc%energy%vdw, ks%calc%vxc)
989 if (ks%calc%calc_energy)
then
1008 subroutine v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
1009 type(
v_ks_t),
target,
intent(inout) :: ks
1012 class(
space_t),
intent(in) :: space
1016 logical,
optional,
intent(in) :: force_semilocal
1018 integer :: ip, ispin
1021 real(real64) :: exx_energy
1022 real(real64) :: factor
1026 assert(ks%calc%calculating)
1027 ks%calc%calculating = .false.
1029 if (ks%frozen_hxc)
then
1035 safe_deallocate_a(hm%energy)
1036 call move_alloc(ks%calc%energy, hm%energy)
1038 if (hm%self_induced_magnetic)
then
1039 hm%a_ind(1:ks%gr%np, 1:space%dim) = ks%calc%a_ind(1:ks%gr%np, 1:space%dim)
1040 hm%b_ind(1:ks%gr%np, 1:space%dim) = ks%calc%b_ind(1:ks%gr%np, 1:space%dim)
1042 safe_deallocate_a(ks%calc%a_ind)
1043 safe_deallocate_a(ks%calc%b_ind)
1046 if (
allocated(hm%v_static))
then
1047 hm%energy%intnvstatic =
dmf_dotp(ks%gr, ks%calc%total_density, hm%v_static)
1049 hm%energy%intnvstatic =
m_zero
1055 hm%energy%intnvxc =
m_zero
1056 hm%energy%hartree =
m_zero
1057 hm%energy%exchange =
m_zero
1058 hm%energy%exchange_hf =
m_zero
1059 hm%energy%correlation =
m_zero
1062 hm%energy%hartree =
m_zero
1063 call v_ks_hartree(namespace, ks, space, hm, ext_partners)
1069 call dxc_oep_calc(ks%sic%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1070 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1072 call zxc_oep_calc(ks%sic%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1073 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1082 call dxc_oep_calc(ks%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1083 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1085 call zxc_oep_calc(ks%oep, namespace, ks%xc, ks%gr, hm, st, space, &
1086 latt%rcell_volume, hm%energy%exchange, hm%energy%correlation, vxc = ks%calc%vxc)
1092 if (ks%calc%calc_energy)
then
1094 hm%energy%intnvxc =
m_zero
1097 do ispin = 1, hm%d%nspin
1098 if (ispin <= 2)
then
1103 hm%energy%intnvxc = hm%energy%intnvxc + &
1104 factor*
dmf_dotp(ks%gr, st%rho(:, ispin), ks%calc%vxc(:, ispin), reduce = .false.)
1106 if (ks%gr%parallel_in_domains)
call ks%gr%allreduce(hm%energy%intnvxc)
1111 if (ks%theory_level /=
hartree .and. ks%theory_level /=
rdmft)
then
1113 safe_deallocate_a(hm%vxc)
1114 call move_alloc(ks%calc%vxc, hm%vxc)
1117 call lalg_copy(ks%gr%np, hm%d%nspin, ks%calc%vtau, hm%vtau)
1118 call hm%hm_base%accel_copy_pot(ks%gr, hm%vtau)
1119 safe_deallocate_a(ks%calc%vtau)
1125 hm%energy%intnvxc = hm%energy%intnvxc &
1128 hm%energy%intnvxc = hm%energy%intnvxc &
1138 if (.not. ks%xc_photon_include_hartree)
then
1139 hm%energy%hartree =
m_zero
1146 hm%vhxc(ip, 1) = hm%vxc(ip, 1) + hm%vhartree(ip)
1148 if (
allocated(hm%vberry))
then
1150 hm%vhxc(ip, 1) = hm%vhxc(ip, 1) + hm%vberry(ip, 1)
1156 hm%vhxc(ip, 2) = hm%vxc(ip, 2) + hm%vhartree(ip)
1158 if (
allocated(hm%vberry))
then
1160 hm%vhxc(ip, 2) = hm%vhxc(ip, 2) + hm%vberry(ip, 2)
1165 if (hm%d%ispin ==
spinors)
then
1168 hm%vhxc(ip, ispin) = hm%vxc(ip, ispin)
1174 hm%energy%exchange_hf =
m_zero
1179 if (
associated(hm%exxop%st))
then
1182 safe_deallocate_p(hm%exxop%st)
1184 if (
associated(ks%calc%hf_st) .and. hm%exxop%useACE)
then
1200 ks%calc%hf_st, xst, hm%kpoints, exx_energy)
1204 ks%calc%hf_st, xst, hm%kpoints, exx_energy)
1205 if (hm%phase%is_allocated())
then
1212 exx_energy = exx_energy + hm%exxop%singul%energy
1216 select case (ks%theory_level)
1220 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1225 hm%energy%exchange_hf = hm%energy%exchange_hf + exx_energy
1244 if (ks%has_photons .and. (ks%xc_photon == 0))
then
1245 if (
associated(ks%pt_mx%vmf))
then
1246 forall(ip = 1:ks%gr%np) hm%vhxc(ip, 1) = hm%vhxc(ip, 1) + ks%pt_mx%vmf(ip)
1248 forall(ip = 1:ks%gr%np) hm%vhxc(ip, 2) = hm%vhxc(ip, 2) + ks%pt_mx%vmf(ip)
1251 hm%ep%photon_forces(1:space%dim) = ks%pt_mx%fmf(1:space%dim)
1254 if (ks%vdw%vdw_correction /= option__vdwcorrection__none)
then
1255 hm%ep%vdw_forces = ks%vdw%forces
1256 hm%ep%vdw_stress = ks%vdw%stress
1257 safe_deallocate_a(ks%vdw%forces)
1259 hm%ep%vdw_forces = 0.0_real64
1262 if (ks%calc%time_present .or. hm%time_zero)
then
1263 call hm%update(ks%gr, namespace, space, ext_partners, time = ks%calc%time)
1269 safe_deallocate_a(ks%calc%density)
1270 if (ks%calc%total_density_alloc)
then
1271 safe_deallocate_p(ks%calc%total_density)
1273 nullify(ks%calc%total_density)
1284 subroutine v_ks_hartree(namespace, ks, space, hm, ext_partners)
1286 type(
v_ks_t),
intent(inout) :: ks
1287 class(
space_t),
intent(in) :: space
1295 call dpoisson_solve(hm%psolver, namespace, hm%vhartree, ks%calc%total_density)
1301 if (ks%calc%calc_energy)
then
1303 hm%energy%hartree =
m_half*
dmf_dotp(ks%gr, ks%calc%total_density, hm%vhartree)
1307 if(ks%calc%time_present)
then
1310 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick, time=ks%calc%time)
1313 ks%calc%total_density, hm%energy%pcm_corr, time=ks%calc%time)
1318 ks%calc%total_density, hm%energy%pcm_corr, kick=hm%kick)
1321 ks%calc%total_density, hm%energy%pcm_corr)
1332 type(
v_ks_t),
intent(inout) :: ks
1336 ks%frozen_hxc = .
true.
1343 type(
v_ks_t),
intent(inout) :: this
1344 logical,
intent(in) :: calc_cur
1348 this%calculate_current = calc_cur
constant times a vector plus a vector
Copies a vector x, to a vector y.
scales a vector by a constant
This is the common interface to a sorting routine. It performs the shell algorithm,...
pure logical function, public accel_is_enabled()
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
subroutine, public current_init(this, namespace)
type(debug_t), save, public debug
This module implements a calculator for the density and defines related functions.
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
real(real64) function, public zenergy_calc_electronic(namespace, hm, der, st, terms)
real(real64) function, public denergy_calc_electronic(namespace, hm, der, st, terms)
subroutine, public energy_calc_eigenvalues(namespace, hm, der, st)
subroutine, public energy_copy(ein, eout)
subroutine, public dexchange_operator_ace(this, namespace, mesh, st, xst, phase)
subroutine, public zexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, ex, F_out)
subroutine, public dexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, ex, F_out)
subroutine, public zexchange_operator_ace(this, namespace, mesh, st, xst, phase)
subroutine, public exchange_operator_reinit(this, omega, alpha, beta, st)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
integer, parameter, public term_mgga
integer, parameter, public term_dft_u
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public rdmft
integer, parameter, public hartree
logical function, public hamiltonian_elec_has_kick(hm)
integer, parameter, public hartree_fock
integer, parameter, public independent_particles
logical function, public hamiltonian_elec_needs_current(hm, states_are_real)
subroutine, public hamiltonian_elec_update_pot(this, mesh, accumulate)
Update the KS potential of the electronic Hamiltonian.
integer, parameter, public kohn_sham_dft
This module defines classes and functions for interaction partners.
integer, parameter, public dft_u_none
subroutine, public magnetic_constrain_update(this, mesh, std, space, latt, pos)
subroutine, public magnetic_induced(namespace, gr, st, psolver, kpoints, a_ind, b_ind)
This subroutine receives as input a current, and produces as an output the vector potential that it i...
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
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 messages_experimental(name, namespace)
This module handles the communicators for the various parallelization strategies.
logical function, public parse_is_defined(namespace, name)
subroutine, public pcm_hartree_potential(pcm, space, mesh, psolver, ext_partners, vhartree, density, pcm_corr, kick, time)
PCM reaction field due to the electronic density.
subroutine, public mf_calc(this, gr, st, ions, pt_mode, time)
subroutine, public dpoisson_solve_start(this, rho)
subroutine, public dpoisson_solve_finish(this, pot)
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel)
Calculates the Poisson equation. Given the density returns the corresponding potential.
logical pure function, public poisson_is_async(this)
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.
integer, parameter, public pseudo_exchange_unknown
integer, parameter, public pseudo_correlation_unknown
integer, parameter, public pseudo_correlation_any
integer, parameter, public pseudo_exchange_any
This module is intended to contain "only mathematical" functions and procedures.
integer, parameter, private libxc_c_index
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
subroutine, public states_elec_allocate_current(st, space, mesh)
This module provides routines for communicating states when using states parallelization.
subroutine, public states_elec_parallel_remote_access_stop(this)
stop remote memory access for states on other processors
subroutine, public states_elec_parallel_remote_access_start(this)
start remote memory access for states on other processors
subroutine v_ks_hartree(namespace, ks, space, hm, ext_partners)
Hartree contribution to the KS potential. This function is designed to be used by v_ks_calc_finish an...
subroutine, public v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
subroutine, public v_ks_freeze_hxc(ks)
subroutine, public v_ks_end(ks)
subroutine, public v_ks_calculate_current(this, calc_cur)
subroutine, public v_ks_write_info(ks, iunit, namespace)
subroutine, public v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, calc_energy, calc_current, force_semilocal)
This routine starts the calculation of the Kohn-Sham potential. The routine v_ks_calc_finish must be ...
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
subroutine, public v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
subroutine, public x_slater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Interface to X(slater_calc)
subroutine, public x_fbe_calc(id, namespace, psolver, gr, st, space, ex, vxc)
Interface to X(x_fbe_calc) Two possible run modes possible: adiabatic and Sturm-Liouville....
subroutine, public fbe_c_lda_sl(namespace, psolver, gr, st, space, ec, vxc)
Sturm-Liouville version of the FBE local-density correlation functional.
integer, parameter, public xc_family_ks_inversion
declaring 'family' constants for 'functionals' not handled by libxc careful not to use a value define...
integer function, public xc_get_default_functional(dim, pseudo_x_functional, pseudo_c_functional)
Returns the default functional given the one parsed from the pseudopotentials and the space dimension...
integer, parameter, public xc_family_nc_mgga
integer, parameter, public xc_oep_x
Exact exchange.
integer, parameter, public xc_lda_c_fbe_sl
LDA correlation based ib the force-balance equation - Sturm-Liouville version.
integer, parameter, public xc_family_nc_lda
integer, parameter, public xc_oep_x_fbe_sl
Exchange approximation based on the force balance equation - Sturn-Liouville version.
integer, parameter, public xc_oep_x_fbe
Exchange approximation based on the force balance equation.
integer, parameter, public xc_oep_x_slater
Slater approximation to the exact exchange.
integer, parameter, public func_c
integer, parameter, public func_x
subroutine, public xc_ks_inversion_end(ks_inv)
subroutine, public xc_ks_inversion_write_info(ks_inversion, iunit, namespace)
subroutine, public xc_ks_inversion_init(ks_inv, namespace, gr, ions, st, xc, mc, space, kpoints)
subroutine, public xc_ks_inversion_calc(ks_inversion, namespace, space, gr, hm, ext_partners, st, vxc, time)
subroutine, public xc_write_info(xcs, iunit, namespace)
subroutine, public xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
subroutine, public xc_end(xcs)
subroutine, public xc_get_vxc(gr, xcs, st, kpoints, psolver, namespace, space, rho, ispin, rcell_volume, vxc, ex, ec, deltaxc, vtau, ex_density, ec_density, stress_xc, force_orbitalfree)
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
subroutine, public xc_get_nc_vxc(gr, xcs, st, kpoints, space, namespace, rho, vxc, ex, ec, vtau, ex_density, ec_density)
This routines is similar to xc_get_vxc but for noncollinear functionals, which are not implemented in...
integer, parameter, public oep_type_mgga
integer, parameter, public oep_level_none
the OEP levels
subroutine, public xc_oep_end(oep)
subroutine, public zxc_oep_calc(oep, namespace, xcs, gr, hm, st, space, rcell_volume, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
subroutine, public dxc_oep_calc(oep, namespace, xcs, gr, hm, st, space, rcell_volume, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
subroutine, public xc_oep_write_info(oep, iunit, namespace)
integer, parameter, public oep_type_exx
The different types of OEP that we can work with.
subroutine, public xc_oep_init(oep, namespace, gr, st, mc, space, oep_type)
subroutine, public zxc_oep_photon_calc(oep, namespace, xcs, gr, hm, st, space, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
subroutine, public dxc_oep_photon_calc(oep, namespace, xcs, gr, hm, st, space, ex, ec, vxc)
This file handles the evaluation of the OEP potential, in the KLI or full OEP as described in S....
This module implements the "photon-free" electron-photon exchange-correlation functional.
integer, parameter, public sic_none
no self-interaction correction
subroutine, public xc_sic_write_info(sic, iunit, namespace)
integer, parameter, public sic_adsic
Averaged density SIC.
subroutine, public xc_sic_init(sic, namespace, gr, st, mc, space)
initialize the SIC object
subroutine, public xc_sic_end(sic)
finalize the SIC and, if needed, the included OEP
integer, parameter, public sic_pz_oep
Perdew-Zunger SIC (OEP way)
integer, parameter, public sic_amaldi
Amaldi correction term.
subroutine, public xc_sic_calc_adsic(sic, namespace, space, gr, st, hm, xc, density, vxc, ex, ec)
Computes the ADSIC potential and energy.
A module that takes care of xc contribution from vdW interactions.
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
The states_elec_t class contains all electronic wave functions.
subroutine get_functional_from_pseudos(x_functional, c_functional)
subroutine v_a_xc(hm, force_semilocal)
subroutine calculate_density()