88 real(real64),
allocatable :: dchi(:,:,:,:)
89 complex(real64),
allocatable :: zchi(:,:,:,:)
96 type(states_elec_t),
public,
pointer :: st => null()
98 type(poisson_t) :: psolver
99 type(singularity_t) :: singul
100 logical :: useACE = .false.
104 type(fourier_space_op_t) :: coulb
108 real(real64),
parameter,
private :: TOL_EXX_WEIGHT = 1.0e-3_real64
113 subroutine ace_init(this, namespace, st)
114 class(ACE_t),
intent(inout) :: this
115 type(namespace_t),
intent(in) :: namespace
116 type(states_elec_t),
intent(in ) :: st
133 if (this%nst > st%nst)
then
137 if (this%nst < st%nst)
then
138 call messages_experimental(
"Adaptively-compressed exchange defined with a subset of states", namespace=namespace)
141 if (this%nst * st%smear%el_per_state < st%qtot)
then
142 write(
message(1),
'(a)')
"ACESize should at least equal the number of occupied states."
154 class(ACE_t),
intent(inout) :: this
159 safe_deallocate_a(this%dchi)
160 safe_deallocate_a(this%zchi)
167 type(exchange_operator_t),
intent(inout) :: this
168 type(namespace_t),
intent(in) :: namespace
169 class(space_t),
intent(in) :: space
170 type(states_elec_t),
intent(in) :: st
171 type(derivatives_t),
intent(in) :: der
172 type(multicomm_t),
intent(in) :: mc
173 type(stencil_t),
intent(in) :: stencil
174 type(kpoints_t),
intent(in) :: kpoints
175 type(xc_cam_t),
intent(in) :: cam
194 if (this%useACE)
call this%ace%init(namespace, st)
198 call poisson_init(this%psolver, namespace, space, der, mc, stencil, st%qtot, &
199 force_serial = .
true., verbose = .false.)
201 call poisson_init(this%psolver, namespace, space, der, mc, stencil, st%qtot, &
202 force_serial = .
true., verbose = .false., force_cmplx = .
true.)
215 if (
present(st))
then
229 if (
associated(this%st) .and. .not. this%useACE)
then
232 safe_deallocate_p(this%st)
253 call batch_scal(mesh%np, this%st%occ(:, hpsib%ik), hpsib)
261#include "exchange_operator_inc.F90"
264#include "complex.F90"
265#include "exchange_operator_inc.F90"
scale a batch by a constant or vector
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dexchange_operator_hartree_apply(this, namespace, mesh, st_d, kpoints, exx_coef, psib, hpsib)
subroutine, public dexchange_operator_ace(this, namespace, mesh, st, xst, phase)
subroutine ace_init(this, namespace, st)
Initialize an instance of ACE_t.
subroutine, public zexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, F_out)
subroutine, public zexchange_operator_commute_r(this, namespace, mesh, st_d, ik, psi, gpsi)
subroutine, public exchange_operator_init(this, namespace, space, st, der, mc, stencil, kpoints, cam)
subroutine, public exchange_operator_reinit(this, cam, st)
subroutine, public dexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, F_out)
subroutine, public zexchange_operator_single(this, namespace, space, mesh, st_d, kpoints, phase, ist, ik, psi, hpsi, rdmft, force_noace)
subroutine, public dexchange_operator_single(this, namespace, space, mesh, st_d, kpoints, phase, ist, ik, psi, hpsi, rdmft, force_noace)
subroutine ace_end(this)
End an instance of ACE_t.
subroutine, public dexchange_operator_commute_r(this, namespace, mesh, st_d, ik, psi, gpsi)
subroutine, public zexchange_operator_hartree_apply(this, namespace, mesh, st_d, kpoints, exx_coef, psib, hpsib)
subroutine, public exchange_operator_end(this)
subroutine, public zexchange_operator_apply(this, namespace, space, mesh, st_d, kpoints, phase, psib, hpsib, rdmft, force_noace)
subroutine, public dexchange_operator_apply(this, namespace, space, mesh, st_d, kpoints, phase, psib, hpsib, rdmft, force_noace)
subroutine, public zexchange_operator_ace(this, namespace, mesh, st, xst, phase)
real(real64) function, public dexchange_operator_compute_ex(mesh, st, xst)
Compute the exact exchange energy.
subroutine, public exchange_operator_rdmft_occ_apply(this, mesh, hpsib)
real(real64) function, public zexchange_operator_compute_ex(mesh, st, xst)
Compute the exact exchange energy.
subroutine, public fourier_space_op_end(this)
This module implements the underlying real-space grid.
This module is intended to contain "only mathematical" functions and procedures.
This module defines functions over batches of mesh functions.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
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_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
This module handles the communicators for the various parallelization strategies.
Some general things and nomenclature:
subroutine, public poisson_init(this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx)
subroutine, public poisson_end(this)
This module is an helper to perform ring-pattern communications among all states.
subroutine, public singularity_end(this)
subroutine, public singularity_init(this, namespace, space, st, kpoints)
pure logical function, public states_are_real(st)
This module provides routines for communicating all batches in a ring-pattern scheme.
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
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
This module defines stencils used in Octopus.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Coulomb-attenuating method parameters, used in the partitioning of the Coulomb potential into a short...