88 complex(real64),
allocatable :: phase(:, :)
91 complex(real64),
public,
allocatable :: phase_corr(:,:)
94 complex(real64),
allocatable :: phase_spiral(:,:)
97 type(accel_mem_t) :: buff_phase
98 type(accel_mem_t) :: buff_phase_spiral
99 type(accel_mem_t),
public :: buff_phase_corr
100 integer :: buff_phase_qn_start
101 real(real64),
public,
pointer :: spin(:,:,:) => null()
129 class(phase_t),
intent(inout) :: phase
130 class(mesh_t),
intent(in) :: gr
131 type(distributed_t),
intent(in) :: kpt
132 type(kpoints_t),
intent(in) :: kpoints
133 type(states_elec_dim_t),
intent(in) :: d
134 type(space_t),
intent(in) :: space
136 integer :: ip, ik, sp
137 integer(int64) :: ip_inner_global
138 real(real64) :: kpoint(space%dim), x_global(space%dim)
145 phase%buff_phase_qn_start = kpt%start
147 if(kpoints%gamma_only())
then
152 safe_allocate(phase%phase(1:gr%np_part, kpt%start:kpt%end))
153 safe_allocate(phase%phase_corr(gr%np+1:gr%np_part, kpt%start:kpt%end))
155 do ik = kpt%start, kpt%end
157 do ip = gr%np + 1, gr%np_part
158 phase%phase_corr(ip, ik) =
m_one
167 if (gr%der%boundaries%spiralBC)
then
169 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
175 safe_allocate(phase%phase_spiral(1:gr%np_part-sp, 1:2))
178 do ip = sp + 1, gr%np_part
182 phase%phase_spiral(ip-sp, 1) = &
183 exp(
m_zi * sum((gr%x(1:space%dim, ip)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
184 phase%phase_spiral(ip-sp, 2) = &
185 exp(-
m_zi * sum((gr%x(1:space%dim, ip)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
190 call accel_write_buffer(phase%buff_phase_spiral, gr%np_part-sp, 2, phase%phase_spiral)
201 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
204 do ik = kpt%start, kpt%end
205 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
207 do ip = 1, gr%np_part
208 phase%phase(ip, ik) =
exp(-
m_zi * sum(gr%x(1:space%dim, ip) * kpoint(1:space%dim)))
214 do ip = sp + 1, gr%np_part
220 phase%phase_corr(ip, ik) = phase%phase(ip, ik)* &
221 exp(
m_zi * sum(x_global(1:space%dim) * kpoint(1:space%dim)))
231 call accel_write_buffer(phase%buff_phase_corr, gr%np_part - gr%np, kpt%nlocal, phase%phase_corr)
239 subroutine phase_update_phases(phase, mesh, kpt, kpoints, d, space, uniform_vector_potential)
240 class(
phase_t),
intent(inout) :: phase
241 class(
mesh_t),
intent(in) :: mesh
245 type(
space_t),
intent(in) :: space
246 real(real64),
allocatable,
intent(in) :: uniform_vector_potential(:)
248 integer :: ik, ip, sp
249 integer(int64),
dimension(2) :: np, gsize, bsize
250 integer(int64) :: ip_inner_global
251 real(real64) :: kpoint(space%dim)
252 real(real64),
allocatable :: x_global(:,:), kpt_vec_pot(:,:)
253 type(
accel_mem_t) :: buff_vec_pot, buff_x_global, buff_x
255 real(real64) :: tmp_sum
257 if (.not.
allocated(uniform_vector_potential))
return
262 if (.not.
allocated(phase%phase))
then
263 safe_allocate(phase%phase(1:mesh%np_part, kpt%start:kpt%end))
266 mesh%np_part*kpt%nlocal)
270 if (.not.
allocated(phase%phase_corr))
then
271 safe_allocate(phase%phase_corr(mesh%np+1:mesh%np_part, kpt%start:kpt%end))
274 (mesh%np_part - mesh%np)*kpt%nlocal)
283 if (mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
285 safe_allocate(x_global(1:space%dim,(sp + 1):mesh%np_part))
288 do ip = sp + 1, mesh%np_part
292 x_global(:,ip) =
mesh_x_global(mesh, ip_inner_global) - mesh%x(1:space%dim, ip)
298 do ik = kpt%start, kpt%end
299 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
301 kpoint(1:space%dim) = kpoint(1:space%dim) + uniform_vector_potential(1:space%dim)
304 do ip = 1, mesh%np_part
305 tmp_sum = sum(mesh%x(1:space%dim, ip)*kpoint(1:space%dim))
306 phase%phase(ip, ik) = cmplx(
cos(tmp_sum), -
sin(tmp_sum), real64)
311 do ip = sp + 1, mesh%np_part
312 tmp_sum = sum(x_global(1:space%dim, ip)*kpoint(1:space%dim))
313 phase%phase_corr(ip, ik) = cmplx(
cos(tmp_sum),
sin(tmp_sum), real64)
322 safe_allocate(kpt_vec_pot(1:space%dim,kpt%start:kpt%end))
323 do ik = kpt%start, kpt%end
324 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
325 kpt_vec_pot(1:space%dim, ik) = kpoint(1:space%dim) + uniform_vector_potential(1:space%dim)
334 call accel_write_buffer(buff_x_global, space%dim, mesh%np_part-sp, x_global(1:space%dim,(sp + 1):mesh%np_part), async=.
true.)
351 np = (/mesh%np_part, kpt%nlocal/)
358 call accel_read_buffer(phase%buff_phase_corr, mesh%np_part - mesh%np, kpt%nlocal, phase%phase_corr)
363 safe_deallocate_a(kpt_vec_pot)
366 safe_deallocate_a(x_global)
374 class(
phase_t),
intent(inout) :: phase
387 safe_deallocate_a(phase%phase)
388 safe_deallocate_a(phase%phase_corr)
389 safe_deallocate_a(phase%phase_spiral)
397 class(
phase_t),
intent(inout) :: phase
398 class(
mesh_t),
intent(in) :: mesh
410 phase%buff_phase_qn_start = kpt%start
416 if (
allocated(phase%phase))
then
417 assert(
size(phase%phase, 1) == mesh%np_part)
418 nlocal = ubound(phase%phase, dim=2) - lbound(phase%phase, dim=2) + 1
423 if (
allocated(phase%phase_corr))
then
424 assert(
size(phase%phase_corr, 1) == mesh%np_part - mesh%np)
425 nlocal = ubound(phase%phase_corr, dim=2) - lbound(phase%phase_corr, dim=2) + 1
427 size(phase%phase_corr, 1)*nlocal)
428 call accel_write_buffer(phase%buff_phase_corr,
size(phase%phase_corr, 1), nlocal, phase%phase_corr)
431 if (
allocated(phase%phase_spiral))
then
433 size(phase%phase_spiral, 1)*
size(phase%phase_spiral, 2))
435 size(phase%phase_spiral, 2), phase%phase_spiral)
445 class(
phase_t),
intent(in) :: phase
446 class(
mesh_t),
intent(in) :: mesh
448 logical,
optional,
intent(in) :: async
451 logical :: phase_correction
456 phase_correction = phase%is_allocated()
460 if (phase_correction)
then
461 call phase%apply_to(mesh, mesh%np, .false., psib, async=async)
471 class(
phase_t),
intent(in) :: phase
472 class(
mesh_t),
intent(in) :: mesh
474 logical,
optional,
intent(in) :: async
476 logical :: phase_correction
481 phase_correction = phase%is_allocated()
485 if (phase_correction)
then
486 call phase%apply_to(mesh, mesh%np, .
true., psib, async=async)
496 class(
phase_t),
intent(in) :: this
497 class(
mesh_t),
intent(in) :: mesh
498 integer,
intent(in) :: np
499 logical,
intent(in) :: conjugate
500 type(
wfs_elec_t),
target,
intent(inout) :: psib
501 type(
wfs_elec_t),
optional,
target,
intent(in) :: src
502 logical,
optional,
intent(in) :: async
504 integer :: ip, ii, sp
506 complex(real64) :: phase
507 integer(int64),
dimension(3) :: gsizes, bsizes
515 assert(np <= mesh%np_part)
517 assert(psib%ik >= lbound(this%phase, dim=2))
518 assert(psib%ik <= ubound(this%phase, dim=2))
521 if (
present(src)) src_ => src
523 assert(src_%has_phase .eqv. conjugate)
524 assert(src_%ik == psib%ik)
528 sp = min(np, mesh%np)
529 if (np > mesh%np .and. mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
531 select case (psib%status())
538 do ip = 1, min(mesh%np, np)
539 phase = conjg(this%phase(ip, psib%ik))
541 do ii = 1, psib%nst_linear
542 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
550 phase = conjg(this%phase(ip, psib%ik))
552 do ii = 1, psib%nst_linear
553 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
562 do ip = 1, min(mesh%np, np)
563 phase = this%phase(ip, psib%ik)
565 do ii = 1, psib%nst_linear
566 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
574 phase = this%phase(ip, psib%ik)
576 do ii = 1, psib%nst_linear
577 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
589 do ii = 1, psib%nst_linear
591 do ip = 1, min(mesh%np, np)
592 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
599 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
607 do ii = 1, psib%nst_linear
609 do ip = 1, min(mesh%np, np)
610 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
617 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
651 psib%has_phase = .not. conjugate
664 class(
phase_t),
intent(in) :: this
665 complex(real64),
intent(inout) :: psi(:, :)
666 integer,
intent(in) :: np
667 integer,
intent(in) :: dim
668 integer,
intent(in) :: ik
669 logical,
intent(in) :: conjugate
675 assert(ik >= lbound(this%phase, dim=2))
676 assert(ik <= ubound(this%phase, dim=2))
685 psi(ip, idim) = conjg(this%phase(ip, ik))*psi(ip, idim)
694 psi(ip, idim) = this%phase(ip, ik)*psi(ip, idim)
710 class(
phase_t),
intent(in) :: this
714 integer :: ip, ii, sp
715 integer,
allocatable :: spin_label(:)
717 integer(int64) :: bsize
718 integer(int64),
dimension(2) :: np, gsizes, bsizes
725 assert(der%boundaries%spiral)
729 if (der%mesh%parallel_in_domains) sp = der%mesh%np + der%mesh%pv%np_ghost
732 select case (psib%status())
736 do ip = sp + 1, der%mesh%np_part
737 do ii = 1, psib%nst_linear, 2
738 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0)
then
739 psib%zff_pack(ii+1, ip) = psib%zff_pack(ii+1, ip)*this%phase_spiral(ip-sp, 1)
741 psib%zff_pack(ii, ip) = psib%zff_pack(ii, ip)*this%phase_spiral(ip-sp, 2)
750 do ii = 1, psib%nst_linear, 2
751 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0)
then
753 do ip = sp + 1, der%mesh%np_part
754 psib%zff_linear(ip, ii+1) = psib%zff_linear(ip, ii+1)*this%phase_spiral(ip-sp, 1)
759 do ip = sp + 1, der%mesh%np_part
760 psib%zff_linear(ip, ii) = psib%zff_linear(ip, ii)*this%phase_spiral(ip-sp, 2)
776 safe_allocate(spin_label(1:psib%nst_linear))
778 do ii = 1, psib%nst_linear, 2
779 if (this%spin(3, psib%linear_to_ist(ii), psib%ik) > 0) spin_label(ii)=1
797 np = (/psib%pack_size(1)/2_int64, int(der%mesh%np_part - sp, int64)/)
798 bsizes = (/psib%pack_size(1)/2, 2*bsize/)
807 safe_deallocate_a(spin_label)
817 logical pure function phase_is_allocated(this)
818 class(
phase_t),
intent(in) :: this
820 phase_is_allocated =
allocated(this%phase)
831 class(
phase_t),
intent(in) :: phase
832 type(grid_t),
intent(in) :: gr
833 type(distributed_t),
intent(in) :: kpt
834 type(wfs_elec_t),
intent(in) :: psib
835 type(wfs_elec_t),
intent(out) :: psib_with_phase
837 integer :: k_offset, n_boundary_points
841 call psib%copy_to(psib_with_phase)
842 if (phase%is_allocated())
then
843 call phase%apply_to(gr, gr%np, conjugate = .false., psib = psib_with_phase, src = psib, async=.
true.)
846 k_offset = psib%ik - kpt%start
847 n_boundary_points = int(gr%np_part - gr%np)
848 call boundaries_set(gr%der%boundaries, gr, psib_with_phase, phase_correction = phase%phase_corr(:, psib%ik), &
849 buff_phase_corr = phase%buff_phase_corr, offset=k_offset*n_boundary_points, async=.
true.)
851 call psib%copy_data_to(gr%np, psib_with_phase)
852 call boundaries_set(gr%der%boundaries, gr, psib_with_phase)
855 call psib_with_phase%do_pack(copy = .
true.)
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
integer function, public accel_kernel_block_size(kernel)
subroutine, public accel_free_buffer(this, async)
subroutine, public accel_kernel_start_call(this, file_name, kernel_name, flags)
subroutine, public accel_finish()
subroutine, public accel_detach_buffer(this)
Clear a buffer handle without freeing device memory.
integer, parameter, public accel_mem_read_write
type(accel_kernel_t), target, save, public kernel_phase_spiral
pure logical function, public accel_is_enabled()
integer, parameter, public accel_mem_read_only
This module implements batches of mesh functions.
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
complex(real64), parameter, public m_zi
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
integer(int64) function, public mesh_periodic_point(mesh, space, ip)
This function returns the point inside the grid corresponding to a boundary point when PBCs are used....
real(real64) function, dimension(1:mesh%box%dim), public mesh_x_global(mesh, ipg)
Given a global point index, this function returns the coordinates of the point.
subroutine phase_phase_spiral(this, der, psib)
apply spiral phase
subroutine phase_unset_phase_corr(phase, mesh, psib, async)
unset the phase correction (if necessary)
subroutine, public phase_accel_rebuild(phase, mesh, kpt)
Rebuild phase accelerator buffers after an intrinsic copy.
subroutine phase_copy_and_set_phase(phase, gr, kpt, psib, psib_with_phase)
Copy a batch to another batch and apply the Bloch phase to it.
subroutine phase_init_phases(phase, gr, kpt, kpoints, d, space)
Initiliaze the phase arrays and copy to GPU the data.
subroutine phase_end(phase)
Releases the memory of the phase object.
subroutine phase_update_phases(phase, mesh, kpt, kpoints, d, space, uniform_vector_potential)
Update the phases.
logical pure function phase_is_allocated(this)
subroutine phase_apply_batch(this, mesh, np, conjugate, psib, src, async)
apply (remove) the phase to the wave functions before (after) applying the Hamiltonian
subroutine phase_set_phase_corr(phase, mesh, psib, async)
set the phase correction (if necessary)
subroutine phase_apply_mf(this, psi, np, dim, ik, conjugate)
apply (or remove) the phase to a wave function psi
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.
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
type(type_t), public type_cmplx
type(type_t), public type_integer
class representing derivatives
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
A container for the phase.
class for organizing spins and k-points
batches of electronic states