34 use,
intrinsic :: iso_fortran_env
70 logical :: needs_vtau = .false.
72 real(real64),
public,
allocatable :: vhartree(:)
73 real(real64),
public,
allocatable :: vxc(:,:)
74 real(real64),
public,
allocatable :: vhxc(:,:)
77 real(real64),
allocatable :: vtau(:,:)
78 real(real64),
allocatable :: lapl_vtau(:,:)
80 type(accel_mem_t) :: vtau_accel
81 type(accel_mem_t) :: lapl_vtau_accel
83 type(derivatives_t),
pointer :: der => null()
118 real(real64),
public,
allocatable :: vhxc(:,:)
119 real(real64),
allocatable :: vtau(:,:)
120 real(real64),
allocatable :: lapl_vtau(:,:)
129 subroutine ks_potential_init(this, der, np, np_part, nspin, theory_level, needs_vtau)
130 class(ks_potential_t),
intent(inout) :: this
131 type(derivatives_t),
target,
intent(in) :: der
132 integer,
intent(in) :: np, np_part
133 integer,
intent(in) :: nspin
134 integer,
intent(in) :: theory_level
135 logical,
intent(in) :: needs_vtau
139 this%theory_level = theory_level
140 this%needs_vtau = needs_vtau
147 safe_allocate(this%vhxc(1:np, 1:nspin))
148 this%vhxc(1:np, 1:nspin) =
m_zero
152 safe_allocate(this%vhartree(1:np_part))
155 safe_allocate(this%vxc(1:np, 1:nspin))
159 safe_allocate(this%vtau(1:np_part, 1:nspin))
161 safe_allocate(this%lapl_vtau(1:np, 1:nspin))
182 safe_deallocate_a(this%vhxc)
183 safe_deallocate_a(this%vhartree)
184 safe_deallocate_a(this%vxc)
185 safe_deallocate_a(this%vtau)
186 safe_deallocate_a(this%lapl_vtau)
199 real(real64),
contiguous,
intent(in) :: vtau(:,:)
203 assert(this%needs_vtau)
216 integer :: offset, ispin
220 assert(this%needs_vtau)
222 do ispin = 1, this%nspin
228 do ispin = 1, this%nspin
229 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
230 call accel_write_buffer(this%lapl_vtau_accel, this%np, this%lapl_vtau(:, ispin), offset = offset)
242 integer :: offset, ispin
254 if (.not. this%needs_vtau)
then
259 assert(
allocated(this%vtau))
260 assert(
allocated(this%lapl_vtau))
268 do ispin = 1, this%nspin
269 call accel_write_buffer(this%vtau_accel, this%np, this%vtau(:, ispin), offset = offset)
270 call accel_write_buffer(this%lapl_vtau_accel, this%np, this%lapl_vtau(:, ispin), offset = offset)
280 real(real64),
contiguous,
intent(inout) :: pot(:,:)
281 integer,
optional,
intent(in) :: nspin
289 assert(
size(pot, dim=1) >= this%np)
290 assert(
size(pot, dim=2) >= nspin_)
304 class(
mesh_t),
intent(in) :: mesh
305 integer,
intent(out) :: ierr
307 integer :: iunit, err, err2(2), isp
308 character(len=MAX_PATH_LEN) :: filename
309 character(len=100) :: lines(2)
321 message(1) =
"Debug: Writing Vhxc restart."
326 iunit = restart%open(
'vhxc')
327 lines(1) =
'# #spin #nspin filename'
329 call restart%write(iunit, lines, 2, err)
330 if (err /= 0) ierr = ierr + 1
333 do isp = 1, this%nspin
335 write(lines(1),
'(i8,a,i8,a)') isp,
' | ', this%nspin,
' | "'//trim(adjustl(filename))//
'"'
336 call restart%write(iunit, lines, 1, err)
337 if (err /= 0) err2(1) = err2(1) + 1
339 call restart%write_mesh_function(filename, mesh, this%vhxc(:,isp), err)
340 if (err /= 0) err2(2) = err2(2) + 1
343 if (err2(1) /= 0) ierr = ierr + 2
344 if (err2(2) /= 0) ierr = ierr + 4
347 call restart%write(iunit, lines, 1, err)
348 if (err /= 0) ierr = ierr + 4
351 if (this%needs_vtau)
then
352 lines(1) =
'# #spin #nspin filename'
354 call restart%write(iunit, lines, 2, err)
355 if (err /= 0) ierr = ierr + 8
358 do isp = 1, this%nspin
360 write(lines(1),
'(i8,a,i8,a)') isp,
' | ', this%nspin,
' | "'//trim(adjustl(filename))//
'"'
361 call restart%write(iunit, lines, 1, err)
362 if (err /= 0) err2(1) = err2(1) + 16
364 call restart%write_mesh_function(filename, mesh, this%vtau(:,isp), err)
365 if (err /= 0) err2(1) = err2(1) + 1
368 if (err2(1) /= 0) ierr = ierr + 32
369 if (err2(2) /= 0) ierr = ierr + 64
372 call restart%write(iunit, lines, 1, err)
373 if (err /= 0) ierr = ierr + 128
376 call restart%close(iunit)
379 message(1) =
"Debug: Writing Vhxc restart done."
392 class(
mesh_t),
intent(in) :: mesh
393 integer,
intent(out) :: ierr
395 integer :: err, err2, isp
396 character(len=MAX_PATH_LEN) :: filename
409 message(1) =
"Debug: Reading Vhxc restart."
414 do isp = 1, this%nspin
417 call restart%read_mesh_function(filename, mesh, this%vhxc(:,isp), err)
418 if (err /= 0) err2 = err2 + 1
421 if (err2 /= 0) ierr = ierr + 1
425 if (this%needs_vtau)
then
426 do isp = 1, this%nspin
429 call restart%read_mesh_function(filename, mesh, this%vtau(:,isp), err)
430 if (err /= 0) err2 = err2 + 1
434 if (err2 /= 0) ierr = ierr + 2
438 message(1) =
"Debug: Reading Vhxc restart done."
449 integer,
optional,
intent(in) :: order
474 real(real64),
intent(in) :: current_time
475 real(real64),
intent(in) :: dt
490 integer,
intent(in) :: history
497 if (.not.
present(storage))
then
498 if (this%needs_vtau)
then
506 if (this%needs_vtau)
then
520 integer,
intent(in) :: history
526 if (this%needs_vtau)
then
540 integer,
intent(in) :: order
541 real(real64),
intent(in) :: current_time
542 real(real64),
intent(in) :: dt
543 real(real64),
intent(in) :: interpolation_time
549 if (this%needs_vtau)
then
551 this%vhxc, vtau = this%vtau)
603 integer(int64),
intent(in) :: how
604 character(len=*),
intent(in) :: dir
605 class(
space_t),
intent(in) :: space
606 class(
mesh_t),
intent(in) :: mesh
607 real(real64),
optional,
intent(in) :: pos(:,:)
608 type(
atom_t),
optional,
intent(in) :: atoms(:)
609 type(
mpi_grp_t),
optional,
intent(in) :: grp
611 character(len=MAX_PATH_LEN) :: fname
617 space, mesh, this%vhartree,
units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
619 do is = 1, this%nspin
622 mesh, this%vxc(:, is),
units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
626 mesh, this%vhxc(:, is),
units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
628 if (this%needs_vtau)
then
631 mesh, this%vtau(:, is),
units_out%energy, err, pos=pos, atoms=atoms, grp = grp)
645 if (.not.
allocated(copy%vhxc))
then
646 safe_allocate(copy%vhxc(1:this%np, 1:this%nspin))
648 if (this%needs_vtau)
then
649 if (.not.
allocated(copy%vtau))
then
650 safe_allocate(copy%vtau(1:this%np, 1:this%nspin))
667 call lalg_copy(this%np, this%nspin, this%vhxc, copy%vhxc)
668 if (this%needs_vtau)
then
669 call lalg_copy(this%np, this%nspin, this%vtau, copy%vtau)
684 assert(
allocated(copy%vhxc))
685 call lalg_copy(this%np, this%nspin, copy%vhxc, this%vhxc)
686 if (this%needs_vtau)
then
687 assert(
allocated(copy%vtau))
688 call this%set_vtau(copy%vtau)
697 integer(int64),
intent(in) :: np
698 integer,
intent(in) :: nspin
699 integer(int64),
intent(in) :: pnp
717 safe_deallocate_a(this%vhxc)
718 safe_deallocate_a(this%vtau)
727 real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot) result(diff)
730 class(
mesh_t),
intent(in) :: mesh
731 real(real64),
intent(in) :: rho(:,:)
732 real(real64),
intent(in) :: qtot
734 real(real64),
allocatable :: diff_pot(:)
737 push_sub(ks_potential_check_convergence)
739 safe_allocate(diff_pot(1:this%np))
744 diff_pot(ip) = abs(copy%vhxc(ip, 1) - this%vhxc(ip, 1))*rho(ip, 1)
746 do is = 2, this%nspin
749 diff_pot(ip) = diff_pot(ip) + abs(copy%vhxc(ip, is) - this%vhxc(ip, is))*rho(ip, is)
755 safe_deallocate_a(diff_pot)
757 pop_sub(ks_potential_check_convergence)
763 type(potential_interpolation_t),
intent(inout) :: vksold
764 real(real64),
intent(in) :: times(:)
765 real(real64),
intent(in) :: current_time
769 assert(
size(times) == 3)
771 call interpolate(times, vksold%v_old(:, :, 1:3), current_time, vksold%v_old(:, :, 0))
772 if (this%needs_vtau)
then
773 call interpolate(times, vksold%vtau_old(:, :, 1:3), current_time, vksold%vtau_old(:, :, 0))
783 real(real64),
intent(in) :: dt
789 assert(not_in_openmp())
792 do ispin = 1, this%nspin
795 vold%vhxc(ip, ispin) = m_half*dt*(this%vhxc(ip, ispin) - vold%vhxc(ip, ispin))
798 if (this%needs_vtau)
then
801 vold%vtau(ip, ispin) = m_half*dt*(this%vtau(ip, ispin) - vold%vtau(ip, ispin))
812#include "complex.F90"
813#include "ks_potential_inc.F90"
817#include "ks_potential_inc.F90"
constant times a vector plus a vector
Copies a vector x, to a vector y.
subroutine, public accel_free_buffer(this, async)
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
pure logical function, public accel_is_enabled()
integer, parameter, public accel_mem_read_only
type(debug_t), save, public debug
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
real(real64), parameter, public m_zero
integer, parameter, public independent_particles
Theory level.
logical pure function, public not_in_openmp()
real(real64), parameter, public m_one
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
A module to handle KS potential, without the external potential.
subroutine ks_potential_set_interpolated_potentials(this, vksold, history)
Set the interpolated potentials to history.
subroutine zks_potential_current_mass_renormalization(this, gpsi, space_dim, ndim, ispin)
Nonlocal contribution of vtau for the current.
subroutine dks_potential_apply_lapl_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of lapl_vtau.
subroutine ks_potential_add_vhxc(this, pot, nspin)
Adds vHxc to the potential.
subroutine dks_potential_current_mass_renormalization(this, gpsi, space_dim, ndim, ispin)
Nonlocal contribution of vtau for the current.
subroutine, public vtau_set_vin(field, this)
subroutine ks_potential_get_interpolated_potentials(this, vksold, history, storage)
Get the interpolated potentials from history.
subroutine zks_potential_apply_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of vtau.
subroutine, public vtau_set_vout(field, this)
subroutine zks_potential_mult_vhxc(this, mf, ispin)
Multiply a mesh function by vHxc.
subroutine ks_potential_init(this, der, np, np_part, nspin, theory_level, needs_vtau)
Allocate the memory for the KS potentials.
subroutine ks_potential_load_vhxc(this, restart, mesh, ierr)
Loads the vhxc potentials.
subroutine ks_potential_interpolate_potentials(this, vksold, order, current_time, dt, interpolation_time)
Interpolate potentials to a new time.
subroutine ks_potential_update_vtau_buffer(this)
Update vtau GPU buffer.
subroutine ks_potential_run_zero_iter(this, vksold)
Run zero iter for the interpolation.
subroutine ks_potential_store_copy(this, copy)
Copy the potentials to a storage object.
subroutine, public ks_potential_accel_rebuild(this)
Rebuild accelerator buffers after an intrinsic copy.
subroutine ks_potential_mix_potentials(this, vold, dt)
Replace vold potentials by 0.5*dt(vold + vhxc)
subroutine ks_potential_perform_interpolation(this, vksold, times, current_time)
Perform a time interpolation of the potentials.
subroutine, public vtau_get_vnew(field, this)
subroutine ks_potential_init_interpolation(this, vksold, order)
Initialize the potential interpolation.
real(real64) function ks_potential_check_convergence(this, copy, mesh, rho, qtot)
Check the convergence of a vhxc for predictor-corrector.
subroutine ks_potential_restore_copy(this, copy)
Copy the potentials from a storage object.
subroutine ks_potential_storage_allocate(this, copy)
Copy the potentials to a storage object.
subroutine ks_potential_dump_vhxc(this, restart, mesh, ierr)
Dumps the vhxc potentials.
subroutine ks_potential_interpolation_new(this, vksold, current_time, dt)
New interpolation point for the interpolation.
subroutine ks_potential_end(this)
Releases the memory for the KS potentials.
subroutine dks_potential_mult_vhxc(this, mf, ispin)
Multiply a mesh function by vHxc.
subroutine, public xc_copied_potentials_end(this)
Finalizer for the copied potentials.
subroutine dks_potential_apply_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of vtau.
subroutine ks_potential_set_vtau(this, vtau)
Set vtau and update the corresponding GPU buffer.
subroutine xc_copied_potentials_copy_vhxc_to_buffer(this, np, nspin, pnp, buffer)
Copy the vhxc potential to a gpu buffer.
subroutine ks_potential_output_potentials(this, namespace, how, dir, space, mesh, pos, atoms, grp)
Outputs vh, vxc, and vtau potentials.
subroutine zks_potential_apply_lapl_vtau_psi(this, mesh, d, ispin, psib, vpsib)
Wrapper to hamiltonian_elec_base_local_sub to hide the data of lapl_vtau.
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
this module contains the low-level part of the output system
character(len=max_path_len) function, public get_filename_with_spin(output, nspin, spin_index)
Returns the filame as output, or output-spX is spin polarized.
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
subroutine, public potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
subroutine, public potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
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.
type(unit_system_t), public units_out
Describes mesh distribution to nodes.
Quantities used in mixing: Input, output and new potentials, and the residuals.
This is defined even when running serial.