66 type(projector_matrix_t),
allocatable,
public :: projector_matrices(:)
67 integer,
public :: nprojector_matrices
68 logical,
public :: apply_projector_matrices
69 logical,
public :: has_non_local_potential
70 integer :: full_projection_size
71 integer,
public :: max_npoints
72 integer,
public :: total_points
74 logical :: projector_mix
75 complex(real64),
allocatable,
public :: projector_phases(:, :, :, :)
76 integer,
allocatable,
public :: projector_to_atom(:)
78 integer,
allocatable :: regions(:)
79 integer,
public :: nphase
84 type(accel_mem_t) :: buff_offsets
85 type(accel_mem_t) :: buff_matrices
86 type(accel_mem_t) :: buff_maps
87 type(accel_mem_t) :: buff_scals
88 type(accel_mem_t) :: buff_position
89 type(accel_mem_t) :: buff_pos
90 type(accel_mem_t) :: buff_invmap
91 type(accel_mem_t),
public :: buff_projector_phases
92 type(accel_mem_t) :: buff_mix
93 logical :: projector_self_overlap
94 real(real64),
pointer,
public :: spin(:,:,:) => null()
131 real(real64),
allocatable :: dprojection(:, :)
132 complex(real64),
allocatable :: zprojection(:, :)
133 type(accel_mem_t) :: buff_projection
134 type(accel_mem_t) :: buff_spin_to_phase
135 type(accel_mem_t),
allocatable :: buff_phasepsi(:)
136 type(accel_mem_t),
allocatable :: buff_projection_temp(:)
137 integer :: cuda_stream_projection_DtH
146 class(nonlocal_pseudopotential_t),
intent(inout) :: this
150 this%apply_projector_matrices = .false.
151 this%has_non_local_potential = .false.
152 this%nprojector_matrices = 0
154 this%projector_self_overlap = .false.
168 integer :: imat, iatom
173 if (.not.
allocated(this%projector_matrices) .or. this%nprojector_matrices == 0)
then
178 assert(
allocated(this%projector_to_atom))
179 assert(
size(this%projector_matrices) >= this%nprojector_matrices)
180 assert(
size(this%projector_to_atom) >= this%nprojector_matrices)
182 do imat = 1, this%nprojector_matrices
183 pmat => this%projector_matrices(imat)
184 iatom = this%projector_to_atom(imat)
186 assert(iatom >= 1 .and. iatom <= epot%natoms)
187 assert(
allocated(epot%proj))
188 assert(
allocated(epot%proj(iatom)%sphere%map))
189 assert(
allocated(epot%proj(iatom)%sphere%rel_x))
191 pmat%map => epot%proj(iatom)%sphere%map
192 pmat%position => epot%proj(iatom)%sphere%rel_x
208 if (
allocated(this%projector_matrices))
then
219 if (
allocated(this%projector_phases))
call accel_free_buffer(this%buff_projector_phases)
222 do iproj = 1, this%nprojector_matrices
225 safe_deallocate_a(this%regions)
226 safe_deallocate_a(this%projector_matrices)
227 safe_deallocate_a(this%projector_phases)
228 safe_deallocate_a(this%projector_to_atom)
243 class(
space_t),
intent(in) :: space
244 class(
mesh_t),
intent(in) :: mesh
245 type(
epot_t),
target,
intent(in) :: epot
247 integer :: iatom, iproj, ll, lmax, lloc, mm, ic, jc
248 integer :: nmat, imat, ip, iorder
249 integer :: nregion, jatom, katom, iregion
250 integer,
allocatable :: order(:), head(:), region_count(:)
251 logical,
allocatable :: atom_counted(:)
265 safe_allocate(order(1:epot%natoms))
266 safe_allocate(head(1:epot%natoms + 1))
267 safe_allocate(region_count(1:epot%natoms))
268 safe_allocate(atom_counted(1:epot%natoms))
270 this%projector_self_overlap = .false.
271 atom_counted = .false.
277 nregion = nregion + 1
278 assert(nregion <= epot%natoms)
280 region_count(nregion) = 0
282 do iatom = 1, epot%natoms
283 if (atom_counted(iatom)) cycle
288 assert(
associated(epot%proj(iatom)%sphere%mesh))
289 do jatom = 1, region_count(nregion)
290 katom = order(head(nregion) + jatom - 1)
292 overlap =
submesh_overlap(epot%proj(iatom)%sphere, epot%proj(katom)%sphere, space)
297 if (.not. overlap)
then
300 region_count(nregion) = region_count(nregion) + 1
301 order(head(nregion) - 1 + region_count(nregion)) = iatom
302 atom_counted(iatom) = .
true.
307 head(nregion + 1) = head(nregion) + region_count(nregion)
309 if (all(atom_counted))
exit
312 safe_deallocate_a(atom_counted)
313 safe_deallocate_a(region_count)
320 do iregion = 1, nregion
321 do iatom = head(iregion), head(iregion + 1) - 1
323 do jatom = head(iregion), iatom - 1
325 assert(.not.
submesh_overlap(epot%proj(order(iatom))%sphere, epot%proj(order(jatom))%sphere, space))
336 this%nprojector_matrices = 0
337 this%apply_projector_matrices = .false.
338 this%has_non_local_potential = .false.
339 this%nregions = nregion
342 do iorder = 1, epot%natoms
343 iatom = order(iorder)
346 this%has_non_local_potential = .
true.
351 do iorder = 1, epot%natoms
352 iatom = order(iorder)
355 this%nprojector_matrices = this%nprojector_matrices + 1
356 this%apply_projector_matrices = .
true.
361 if (mesh%use_curvilinear) this%apply_projector_matrices = .false.
363 if (.not. this%apply_projector_matrices)
then
364 safe_deallocate_a(order)
365 safe_deallocate_a(head)
372 safe_allocate(this%projector_matrices(1:this%nprojector_matrices))
373 safe_allocate(this%regions(1:this%nprojector_matrices + 1))
374 safe_allocate(this%projector_to_atom(1:epot%natoms))
376 this%full_projection_size = 0
377 this%regions(this%nregions + 1) = this%nprojector_matrices + 1
379 this%projector_mix = .false.
382 do iregion = 1, this%nregions
383 this%regions(iregion) = iproj + 1
384 do iorder = head(iregion), head(iregion + 1) - 1
386 iatom = order(iorder)
392 pmat => this%projector_matrices(iproj)
394 this%projector_to_atom(iproj) = iatom
396 lmax = epot%proj(iatom)%lmax
397 lloc = epot%proj(iatom)%lloc
404 if (ll == lloc) cycle
406 nmat = nmat + epot%proj(iatom)%kb_p(ll, mm)%n_c
416 if (ll == lloc) cycle
418 kb_p => epot%proj(iatom)%kb_p(ll, mm)
420 call lalg_copy(pmat%npoints, kb_p%p(:, ic), pmat%dprojectors(:, imat))
421 pmat%scal(imat) = kb_p%e(ic)*mesh%vol_pp(1)
427 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
431 this%projector_mix = .
true.
436 if (ll == lloc) cycle
456 if (ll == lloc) cycle
458 hgh_p => epot%proj(iatom)%hgh_p(ll, mm)
465 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 1) = hgh_p%h(ic, jc) +
m_half*mm*hgh_p%k(ic, jc)
466 pmat%zmix(imat - 1 + ic, imat - 1 + jc, 2) = hgh_p%h(ic, jc) -
m_half*mm*hgh_p%k(ic, jc)
469 pmat%zmix(imat - 1 + ic, imat + 3 - 1 + jc, 3) =
m_half*hgh_p%k(ic, jc) * &
470 sqrt(real(ll*(ll+1)-mm*(mm+1), real64))
474 pmat%zmix(imat - 1 + ic, imat - 3 - 1 + jc, 4) =
m_half*hgh_p%k(ic, jc) * &
475 sqrt(real(ll*(ll+1)-mm*(mm-1), real64))
482 pmat%dmix(imat - 1 + ic, imat - 1 + jc) = hgh_p%h(ic, jc)
489 call lalg_copy(pmat%npoints, hgh_p%zp(:, ic), pmat%zprojectors(:, imat))
491 call lalg_copy(pmat%npoints, hgh_p%dp(:, ic), pmat%dprojectors(:, imat))
493 pmat%scal(imat) = mesh%volume_element
500 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
505 this%projector_mix = .
true.
509 if (lloc /= 0) nmat = nmat + epot%proj(iatom)%kb_p(1, 1)%n_c
512 if (ll == lloc) cycle
514 nmat = nmat + epot%proj(iatom)%rkb_p(ll, mm)%n_c
519 has_mix_matrix = .
true., is_cmplx = .
true.)
526 kb_p => epot%proj(iatom)%kb_p(1, 1)
529 pmat%zmix(ic, ic, 1:2) = kb_p%e(ic)
530 do ip = 1, pmat%npoints
531 pmat%zprojectors(ip, ic) = kb_p%p(ip, ic)
533 pmat%scal(ic) = mesh%volume_element
540 if (ll == lloc) cycle
542 rkb_p => epot%proj(iatom)%rkb_p(ll, mm)
545 do ic = 0, rkb_p%n_c/2-1
546 pmat%zmix(imat + ic*2, imat + ic*2, 1) = rkb_p%f(ic*2+1, 1, 1)
547 pmat%zmix(imat + ic*2, imat + ic*2, 2) = rkb_p%f(ic*2+1, 2, 2)
549 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 1) = rkb_p%f(ic*2+2, 1, 1)
550 pmat%zmix(imat + ic*2+1, imat + ic*2+1, 2) = rkb_p%f(ic*2+2, 2, 2)
553 pmat%zmix(imat + ic*2+rkb_p%n_c, imat + ic*2, 4) = rkb_p%f(ic*2+1, 2, 1)
554 pmat%zmix(imat + ic*2+1+rkb_p%n_c, imat + ic*2+1, 4) = rkb_p%f(ic*2+2, 2, 1)
558 pmat%zmix(imat + ic*2-rkb_p%n_c, imat + ic*2, 3) = rkb_p%f(ic*2+1, 1, 2)
559 pmat%zmix(imat + ic*2+1-rkb_p%n_c, imat + ic*2+1, 3) = rkb_p%f(ic*2+2, 1, 2)
564 call lalg_copy(pmat%npoints, rkb_p%ket(:, ic, 1, 1), pmat%zprojectors(:, imat))
565 pmat%scal(imat) = mesh%volume_element
573 this%projector_self_overlap = this%projector_self_overlap .or. epot%proj(iatom)%sphere%overlap
579 pmat%map => epot%proj(iatom)%sphere%map
580 pmat%position => epot%proj(iatom)%sphere%rel_x
582 pmat%regions = epot%proj(iatom)%sphere%regions
584 this%full_projection_size = this%full_projection_size + pmat%nprojs
589 if (mesh%parallel_in_domains)
then
590 call mesh%mpi_grp%allreduce_inplace(this%projector_self_overlap, 1, mpi_logical, mpi_lor)
593 safe_deallocate_a(order)
594 safe_deallocate_a(head)
596 this%total_points = 0
599 do imat = 1, this%nprojector_matrices
600 pmat => this%projector_matrices(imat)
602 this%max_npoints = max(this%max_npoints, pmat%npoints)
603 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
604 this%total_points = this%total_points + pmat%npoints
620 class(
space_t),
intent(in) :: space
621 class(
mesh_t),
intent(in) :: mesh
632 if (.not.
allocated(this%projector_matrices) .or. this%nprojector_matrices <= 0)
then
667 class(
space_t),
intent(in) :: space
668 class(
mesh_t),
intent(in) :: mesh
670 integer,
parameter :: OFFSET_SIZE = 6
671 integer,
parameter :: POINTS = 1, projs = 2, matrix = 3, map = 4, scal = 5, mix = 6
672 integer :: imat, matrix_size, scal_size
673 integer :: ip, is, ii, ipos, mix_offset
674 integer,
allocatable :: cnt(:), invmap(:, :), invmap2(:), pos(:)
675 integer,
allocatable :: offsets(:, :)
680 assert(
allocated(this%projector_matrices))
681 assert(this%nprojector_matrices > 0)
683 safe_allocate(offsets(1:offset_size, 1:this%nprojector_matrices))
684 safe_allocate(cnt(1:mesh%np))
689 this%total_points = 0
694 do imat = 1, this%nprojector_matrices
695 pmat => this%projector_matrices(imat)
697 this%max_npoints = max(this%max_npoints, pmat%npoints)
698 this%max_nprojs = max(this%max_nprojs, pmat%nprojs)
700 offsets(points, imat) = pmat%npoints
701 offsets(projs, imat) = pmat%nprojs
703 offsets(matrix, imat) = matrix_size
704 matrix_size = matrix_size + pmat%npoints*pmat%nprojs
706 offsets(map, imat) = this%total_points
707 this%total_points = this%total_points + pmat%npoints
709 offsets(scal, imat) = scal_size
710 scal_size = scal_size + pmat%nprojs
712 offsets(mix, imat) = mix_offset
713 if (
allocated(pmat%dmix))
then
714 mix_offset = mix_offset + pmat%nprojs**2
715 else if (
allocated(pmat%zmix))
then
716 mix_offset = mix_offset + 4*pmat%nprojs**2
718 offsets(mix, imat) = -1
721 do is = 1, pmat%npoints
723 cnt(ip) = cnt(ip) + 1
727 safe_allocate(invmap(1:max(maxval(cnt), 1), 1:mesh%np))
728 safe_allocate(invmap2(1:max(maxval(cnt)*mesh%np, 1)))
729 safe_allocate(pos(1:mesh%np + 1))
733 do imat = 1, this%nprojector_matrices
734 pmat => this%projector_matrices(imat)
735 do is = 1, pmat%npoints
737 cnt(ip) = cnt(ip) + 1
738 invmap(cnt(ip), ip) = ii
748 invmap2(ipos) = invmap(ii, ip)
753 if (this%projector_matrices(1)%is_cmplx)
then
762 if (mix_offset > 0)
then
763 if (
allocated(this%projector_matrices(1)%zmix))
then
770 do imat = 1, this%nprojector_matrices
771 pmat => this%projector_matrices(imat)
772 if (pmat%npoints > 0)
then
773 if (pmat%is_cmplx)
then
774 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%zprojectors, &
775 offset = offsets(matrix, imat))
777 call accel_write_buffer(this%buff_matrices, pmat%npoints, pmat%nprojs, pmat%dprojectors, &
778 offset = offsets(matrix, imat))
780 call accel_write_buffer(this%buff_maps, pmat%npoints, pmat%map, offset = offsets(map, imat))
782 offset = 3*offsets(map, imat))
784 call accel_write_buffer(this%buff_scals, pmat%nprojs, pmat%scal, offset = offsets(scal, imat))
785 if (offsets(mix, imat) /= -1)
then
786 if (
allocated(pmat%zmix))
then
787 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, 4, pmat%zmix, offset = offsets(mix, imat))
789 call accel_write_buffer(this%buff_mix, pmat%nprojs, pmat%nprojs, pmat%dmix, offset = offsets(mix, imat))
795 call accel_write_buffer(this%buff_offsets, offset_size, this%nprojector_matrices, offsets)
803 safe_deallocate_a(offsets)
804 safe_deallocate_a(cnt)
805 safe_deallocate_a(invmap)
806 safe_deallocate_a(invmap2)
807 safe_deallocate_a(pos)
817 integer :: ik, imat, iphase, nphase, offset, npoints
821 if (.not.
allocated(this%projector_phases))
then
826 nphase =
size(this%projector_phases, 2)
830 this%total_points*nphase*
size(this%projector_phases, 4))
833 do ik = lbound(this%projector_phases, 4), ubound(this%projector_phases, 4)
834 do imat = 1, this%nprojector_matrices
835 npoints = this%projector_matrices(imat)%npoints
836 do iphase = 1, nphase
837 if (npoints > 0)
then
838 call accel_write_buffer(this%buff_projector_phases, npoints, this%projector_phases(1:, iphase, imat, ik), &
839 offset = offset, async=.
true.)
841 offset = offset + npoints
853 logical pure function nonlocal_pseudopotential_self_overlap(this) result(projector_self_overlap)
856 projector_self_overlap = this%projector_self_overlap
861#include "nonlocal_pseudopotential_inc.F90"
864#include "complex.F90"
865#include "nonlocal_pseudopotential_inc.F90"
Copies a vector x, to a vector y.
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public accel_free_buffer(this, async)
subroutine, public accel_finish()
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
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
This module contains interfaces for BLAS routines You should not use these routines directly....
integer, parameter, public spin_orbit
integer, parameter, public fully_relativistic_zora
real(real64), parameter, public m_zero
real(real64), parameter, public m_half
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
subroutine nonlocal_pseudopotential_destroy_proj(this)
Destroy the data of nonlocal_pseudopotential_t.
subroutine nonlocal_pseudopotential_detach_accel_buffers(this)
Clear copied accelerator handles so they do not alias source buffers.
subroutine dnonlocal_pseudopotential_force(this, mesh, st, spiral_bnd, iqn, ndim, psi1b, psi2b, force)
calculate contribution to forces, from non-local potentials
subroutine znonlocal_pseudopotential_position_commutator(this, mesh, std, spiral_bnd, psib, commpsib, async)
apply the commutator between the non-local potential and the position to the wave functions.
subroutine znonlocal_pseudopotential_force(this, mesh, st, spiral_bnd, iqn, ndim, psi1b, psi2b, force)
calculate contribution to forces, from non-local potentials
subroutine dnonlocal_pseudopotential_start(this, mesh, std, spiral_bnd, psib, projection, async)
Start application of non-local potentials (stored in the Hamiltonian) to the wave functions.
subroutine dnonlocal_pseudopotential_finish(this, mesh, spiral_bnd, std, projection, vpsib)
finish the application of non-local potentials.
subroutine, public nonlocal_pseudopotential_accel_rebuild(this, space, mesh)
Rebuild accelerator buffers after an intrinsic copy.
subroutine dnonlocal_pseudopotential_position_commutator(this, mesh, std, spiral_bnd, psib, commpsib, async)
apply the commutator between the non-local potential and the position to the wave functions.
subroutine, public nonlocal_pseudopotential_rebind_projectors(this, epot)
Rebind projector matrix pointers (map, position) to a target epot.
subroutine nonlocal_pseudopotential_build_accel_buffers(this, space, mesh)
Build accelerator buffers for projectors from host-side projector matrices.
subroutine znonlocal_pseudopotential_r_vnlocal(this, mesh, std, spiral_bnd, psib, commpsib)
Accumulates to commpsib the result of x V_{nl} | psib >
logical pure function nonlocal_pseudopotential_self_overlap(this)
Returns .true. if the Hamiltonian contains projectors, which overlap with themself.
subroutine nonlocal_pseudopotential_init(this)
initialize the nonlocal_pseudopotential_t object
subroutine dnonlocal_pseudopotential_r_vnlocal(this, mesh, std, spiral_bnd, psib, commpsib)
Accumulates to commpsib the result of x V_{nl} | psib >
subroutine znonlocal_pseudopotential_start(this, mesh, std, spiral_bnd, psib, projection, async)
Start application of non-local potentials (stored in the Hamiltonian) to the wave functions.
subroutine nonlocal_pseudopotential_build_projector_phase_accel_buffer(this)
Rebuild projector phase accelerator buffer from host-side projector phases.
subroutine nonlocal_pseudopotential_build_proj(this, space, mesh, epot)
build the projectors for the application of pseudo-potentials
subroutine znonlocal_pseudopotential_finish(this, mesh, spiral_bnd, std, projection, vpsib)
finish the application of non-local potentials.
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.
subroutine, public projector_matrix_deallocate(this)
subroutine, public projector_matrix_allocate(this, nprojs, sphere, has_mix_matrix, is_cmplx)
logical elemental function, public projector_is(p, type)
logical elemental function, public projector_is_null(p)
integer, parameter, public proj_hgh
integer, parameter, public proj_rkb
integer, parameter, public proj_none
integer, parameter, public proj_kb
This module handles spin dimensions of the states and the k-point distribution.
logical function, public submesh_overlap(sm1, sm2, space)
type(type_t), public type_float
type(type_t), public type_cmplx
type(type_t), public type_integer
Describes mesh distribution to nodes.
nonlocal part of the pseudopotential
Class for projections of wave functions.
A set of projectors defined on a submesh.
The rkb_projector data type holds the KB projectors build with total angular momentum eigenfunctions....