87 complex(real64),
allocatable :: phase(:, :)
90 complex(real64),
public,
allocatable :: phase_corr(:,:)
93 complex(real64),
allocatable :: phase_spiral(:,:)
96 type(accel_mem_t) :: buff_phase
97 type(accel_mem_t) :: buff_phase_spiral
98 type(accel_mem_t),
public :: buff_phase_corr
99 integer :: buff_phase_qn_start
100 real(real64),
public,
pointer :: spin(:,:,:) => null()
128 class(phase_t),
intent(inout) :: phase
129 class(mesh_t),
intent(in) :: gr
130 type(distributed_t),
intent(in) :: kpt
131 type(kpoints_t),
intent(in) :: kpoints
132 type(states_elec_dim_t),
intent(in) :: d
133 type(space_t),
intent(in) :: space
135 integer :: ip, ik, sp
136 integer(int64) :: ip_inner_global
137 real(real64) :: kpoint(space%dim), x_global(space%dim)
144 phase%buff_phase_qn_start = kpt%start
146 if(kpoints%gamma_only())
then
151 safe_allocate(phase%phase(1:gr%np_part, kpt%start:kpt%end))
152 safe_allocate(phase%phase_corr(gr%np+1:gr%np_part, kpt%start:kpt%end))
154 do ik = kpt%start, kpt%end
156 do ip = gr%np + 1, gr%np_part
157 phase%phase_corr(ip, ik) =
m_one
166 if (gr%der%boundaries%spiralBC)
then
168 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
174 safe_allocate(phase%phase_spiral(1:gr%np_part-sp, 1:2))
177 do ip = sp + 1, gr%np_part
181 phase%phase_spiral(ip-sp, 1) = &
182 exp(
m_zi * sum((gr%x(1:space%dim, ip)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
183 phase%phase_spiral(ip-sp, 2) = &
184 exp(-
m_zi * sum((gr%x(1:space%dim, ip)-x_global(1:space%dim)) * gr%der%boundaries%spiral_q(1:space%dim)))
189 call accel_write_buffer(phase%buff_phase_spiral, gr%np_part-sp, 2, phase%phase_spiral)
200 if (gr%parallel_in_domains) sp = gr%np + gr%pv%np_ghost
203 do ik = kpt%start, kpt%end
204 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
206 do ip = 1, gr%np_part
207 phase%phase(ip, ik) =
exp(-
m_zi * sum(gr%x(1:space%dim, ip) * kpoint(1:space%dim)))
213 do ip = sp + 1, gr%np_part
219 phase%phase_corr(ip, ik) = phase%phase(ip, ik)* &
220 exp(
m_zi * sum(x_global(1:space%dim) * kpoint(1:space%dim)))
230 call accel_write_buffer(phase%buff_phase_corr, gr%np_part - gr%np, kpt%nlocal, phase%phase_corr)
238 subroutine phase_update_phases(phase, mesh, kpt, kpoints, d, space, uniform_vector_potential)
239 class(
phase_t),
intent(inout) :: phase
240 class(
mesh_t),
intent(in) :: mesh
244 type(
space_t),
intent(in) :: space
245 real(real64),
allocatable,
intent(in) :: uniform_vector_potential(:)
247 integer :: ik, ip, sp
248 integer(int64),
dimension(2) :: np, gsize, bsize
249 integer(int64) :: ip_inner_global
250 real(real64) :: kpoint(space%dim)
251 real(real64),
allocatable :: x_global(:,:), kpt_vec_pot(:,:)
252 type(
accel_mem_t) :: buff_vec_pot, buff_x_global, buff_x
254 real(real64) :: tmp_sum
256 if (.not.
allocated(uniform_vector_potential))
return
261 if (.not.
allocated(phase%phase))
then
262 safe_allocate(phase%phase(1:mesh%np_part, kpt%start:kpt%end))
265 mesh%np_part*kpt%nlocal)
269 if (.not.
allocated(phase%phase_corr))
then
270 safe_allocate(phase%phase_corr(mesh%np+1:mesh%np_part, kpt%start:kpt%end))
273 (mesh%np_part - mesh%np)*kpt%nlocal)
282 if (mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
284 safe_allocate(x_global(1:space%dim,(sp + 1):mesh%np_part))
287 do ip = sp + 1, mesh%np_part
291 x_global(:,ip) =
mesh_x_global(mesh, ip_inner_global) - mesh%x(1:space%dim, ip)
297 do ik = kpt%start, kpt%end
298 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
300 kpoint(1:space%dim) = kpoint(1:space%dim) + uniform_vector_potential(1:space%dim)
303 do ip = 1, mesh%np_part
304 tmp_sum = sum(mesh%x(1:space%dim, ip)*kpoint(1:space%dim))
305 phase%phase(ip, ik) = cmplx(
cos(tmp_sum), -
sin(tmp_sum), real64)
310 do ip = sp + 1, mesh%np_part
311 tmp_sum = sum(x_global(1:space%dim, ip)*kpoint(1:space%dim))
312 phase%phase_corr(ip, ik) = cmplx(
cos(tmp_sum),
sin(tmp_sum), real64)
321 safe_allocate(kpt_vec_pot(1:space%dim,kpt%start:kpt%end))
322 do ik = kpt%start, kpt%end
323 kpoint(1:space%dim) = kpoints%get_point(d%get_kpoint_index(ik))
324 kpt_vec_pot(1:space%dim, ik) = kpoint(1:space%dim) + uniform_vector_potential(1:space%dim)
333 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.)
350 np = (/mesh%np_part, kpt%nlocal/)
357 call accel_read_buffer(phase%buff_phase_corr, mesh%np_part - mesh%np, kpt%nlocal, phase%phase_corr)
362 safe_deallocate_a(kpt_vec_pot)
365 safe_deallocate_a(x_global)
373 class(
phase_t),
intent(inout) :: phase
386 safe_deallocate_a(phase%phase)
387 safe_deallocate_a(phase%phase_corr)
388 safe_deallocate_a(phase%phase_spiral)
397 class(
phase_t),
intent(in) :: phase
398 class(
mesh_t),
intent(in) :: mesh
400 logical,
optional,
intent(in) :: async
403 logical :: phase_correction
408 phase_correction = phase%is_allocated()
412 if (phase_correction)
then
413 call phase%apply_to(mesh, mesh%np, .false., psib, async=async)
423 class(
phase_t),
intent(in) :: phase
424 class(
mesh_t),
intent(in) :: mesh
426 logical,
optional,
intent(in) :: async
428 logical :: phase_correction
433 phase_correction = phase%is_allocated()
437 if (phase_correction)
then
438 call phase%apply_to(mesh, mesh%np, .
true., psib, async=async)
448 class(
phase_t),
intent(in) :: this
449 class(
mesh_t),
intent(in) :: mesh
450 integer,
intent(in) :: np
451 logical,
intent(in) :: conjugate
452 type(
wfs_elec_t),
target,
intent(inout) :: psib
453 type(
wfs_elec_t),
optional,
target,
intent(in) :: src
454 logical,
optional,
intent(in) :: async
456 integer :: ip, ii, sp
458 complex(real64) :: phase
459 integer(int64),
dimension(3) :: gsizes, bsizes
467 assert(np <= mesh%np_part)
469 assert(psib%ik >= lbound(this%phase, dim=2))
470 assert(psib%ik <= ubound(this%phase, dim=2))
473 if (
present(src)) src_ => src
475 assert(src_%has_phase .eqv. conjugate)
476 assert(src_%ik == psib%ik)
480 sp = min(np, mesh%np)
481 if (np > mesh%np .and. mesh%parallel_in_domains) sp = mesh%np + mesh%pv%np_ghost
483 select case (psib%status())
490 do ip = 1, min(mesh%np, np)
491 phase = conjg(this%phase(ip, psib%ik))
493 do ii = 1, psib%nst_linear
494 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
502 phase = conjg(this%phase(ip, psib%ik))
504 do ii = 1, psib%nst_linear
505 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
514 do ip = 1, min(mesh%np, np)
515 phase = this%phase(ip, psib%ik)
517 do ii = 1, psib%nst_linear
518 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
526 phase = this%phase(ip, psib%ik)
528 do ii = 1, psib%nst_linear
529 psib%zff_pack(ii, ip) = phase*src_%zff_pack(ii, ip)
541 do ii = 1, psib%nst_linear
543 do ip = 1, min(mesh%np, np)
544 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
551 psib%zff_linear(ip, ii) = conjg(this%phase(ip, psib%ik))*src_%zff_linear(ip, ii)
559 do ii = 1, psib%nst_linear
561 do ip = 1, min(mesh%np, np)
562 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
569 psib%zff_linear(ip, ii) = this%phase(ip, psib%ik)*src_%zff_linear(ip, ii)
603 psib%has_phase = .not. conjugate
616 class(
phase_t),
intent(in) :: this
617 complex(real64),
intent(inout) :: psi(:, :)
618 integer,
intent(in) :: np
619 integer,
intent(in) :: dim
620 integer,
intent(in) :: ik
621 logical,
intent(in) :: conjugate
627 assert(ik >= lbound(this%phase, dim=2))
628 assert(ik <= ubound(this%phase, dim=2))
637 psi(ip, idim) = conjg(this%phase(ip, ik))*psi(ip, idim)
646 psi(ip, idim) = this%phase(ip, ik)*psi(ip, idim)
662 class(
phase_t),
intent(in) :: this
666 integer :: ip, ii, sp
667 integer,
allocatable :: spin_label(:)
669 integer(int64) :: bsize
670 integer(int64),
dimension(2) :: np, gsizes, bsizes
677 assert(der%boundaries%spiral)
681 if (der%mesh%parallel_in_domains) sp = der%mesh%np + der%mesh%pv%np_ghost
684 select case (psib%status())
688 do ip = sp + 1, der%mesh%np_part
689 do ii = 1, psib%nst_linear, 2
690 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0)
then
691 psib%zff_pack(ii+1, ip) = psib%zff_pack(ii+1, ip)*this%phase_spiral(ip-sp, 1)
693 psib%zff_pack(ii, ip) = psib%zff_pack(ii, ip)*this%phase_spiral(ip-sp, 2)
702 do ii = 1, psib%nst_linear, 2
703 if (this%spin(3,psib%linear_to_ist(ii), psib%ik)>0)
then
705 do ip = sp + 1, der%mesh%np_part
706 psib%zff_linear(ip, ii+1) = psib%zff_linear(ip, ii+1)*this%phase_spiral(ip-sp, 1)
711 do ip = sp + 1, der%mesh%np_part
712 psib%zff_linear(ip, ii) = psib%zff_linear(ip, ii)*this%phase_spiral(ip-sp, 2)
728 safe_allocate(spin_label(1:psib%nst_linear))
730 do ii = 1, psib%nst_linear, 2
731 if (this%spin(3, psib%linear_to_ist(ii), psib%ik) > 0) spin_label(ii)=1
749 np = (/psib%pack_size(1)/2_int64, int(der%mesh%np_part - sp, int64)/)
750 bsizes = (/psib%pack_size(1)/2, 2*bsize/)
759 safe_deallocate_a(spin_label)
769 logical pure function phase_is_allocated(this)
770 class(
phase_t),
intent(in) :: this
772 phase_is_allocated =
allocated(this%phase)
783 class(
phase_t),
intent(in) :: phase
784 type(grid_t),
intent(in) :: gr
785 type(distributed_t),
intent(in) :: kpt
786 type(wfs_elec_t),
intent(in) :: psib
787 type(wfs_elec_t),
intent(out) :: psib_with_phase
789 integer :: k_offset, n_boundary_points
793 call psib%copy_to(psib_with_phase)
794 if (phase%is_allocated())
then
795 call phase%apply_to(gr, gr%np, conjugate = .false., psib = psib_with_phase, src = psib, async=.
true.)
798 k_offset = psib%ik - kpt%start
799 n_boundary_points = int(gr%np_part - gr%np)
800 call boundaries_set(gr%der%boundaries, gr, psib_with_phase, phase_correction = phase%phase_corr(:, psib%ik), &
801 buff_phase_corr = phase%buff_phase_corr, offset=k_offset*n_boundary_points, async=.
true.)
803 call psib%copy_data_to(gr%np, psib_with_phase)
804 call boundaries_set(gr%der%boundaries, gr, psib_with_phase)
807 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()
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 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