72 type(submesh_t) :: sphere
73 complex(real64),
allocatable :: phase(:,:)
76 real(real64) :: Ubar, Jbar
80 real(real64),
allocatable :: V_ij(:,:)
81 real(real64),
allocatable :: coulomb_IIJJ(:,:,:,:,:)
82 complex(real64),
allocatable :: zcoulomb_IIJJ(:,:,:,:,:,:,:)
83 integer,
allocatable:: map_os(:)
84 complex(real64),
allocatable :: phase_shift(:,:)
86 real(real64) :: radius
87 class(species_t),
pointer :: spec => null()
90 real(real64),
allocatable :: dorb(:,:,:)
91 complex(real64),
allocatable :: zorb(:,:,:)
92 complex(real64),
allocatable :: eorb_submesh(:,:,:,:)
93 complex(real64),
allocatable :: eorb_mesh(:,:,:,:)
94 integer :: ldorbs, ldorbs_eorb
95 type(accel_mem_t) :: dbuff_orb, zbuff_orb
96 type(accel_mem_t),
allocatable :: buff_eorb (:)
98 logical :: use_submesh
99 logical :: allocated_on_mesh
103 type(poisson_t) :: poisson
111 type(orbitalset_t),
intent(inout) :: this
137 type(orbitalset_t),
intent(inout) :: this
143 safe_deallocate_a(this%phase)
144 safe_deallocate_a(this%dorb)
145 safe_deallocate_a(this%zorb)
146 safe_deallocate_a(this%eorb_submesh)
147 safe_deallocate_a(this%eorb_mesh)
151 safe_deallocate_a(this%V_ij)
152 safe_deallocate_a(this%coulomb_IIJJ)
153 safe_deallocate_a(this%map_os)
154 safe_deallocate_a(this%phase_shift)
159 if(
allocated(this%buff_eorb))
then
160 do ik = lbound(this%buff_eorb, dim=1), ubound(this%buff_eorb, dim=1)
163 safe_deallocate_a(this%buff_eorb)
172 real(real64),
intent(in) :: jj
173 integer,
intent(in) :: ll, nn
188 integer,
intent(in) :: dim
191 logical,
intent(in) :: spin_polarized
192 real(real64),
optional,
allocatable,
intent(in) :: vec_pot(:)
193 real(real64),
optional,
allocatable,
intent(in) :: vec_pot_var(:, :)
194 integer,
optional,
intent(in) :: kpt_max
196 integer :: ns, iq, is, ikpoint, im, idim, kpt_end
197 real(real64) :: kr, kpoint(1:dim), dx(1:dim)
205 if (
present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
207 do iq = kpt%start, kpt_end
209 if (spin_polarized)
then
210 ikpoint = 1 + (iq - 1)/2
218 kpoint(1:dim) = kpoints%get_point(ikpoint)
224 dx = os%sphere%rel_x(1:dim, is) - os%sphere%mesh%x(os%sphere%map(is), 1:dim) + os%sphere%center(1:dim)
225 kr = dot_product(kpoint, dx)
226 if (
present(vec_pot))
then
227 if (
allocated(vec_pot)) kr = kr + dot_product(vec_pot, dx)
230 if (
present(vec_pot_var))
then
231 if (
allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, os%sphere%map(is)) &
232 *(os%sphere%rel_x(1:dim, is)+os%sphere%center))
235 os%phase(is, iq) =
exp(
m_zi*kr)
238 if (.not. os%use_submesh)
then
242 os%eorb_mesh(:, im, idim, iq) =
m_z0
244 os%eorb_mesh(os%sphere%map(is),im,idim,iq) = os%eorb_mesh(os%sphere%map(is),im,idim,iq) &
245 + os%zorb(is, idim, im) * os%phase(is, iq)
255 os%eorb_submesh(is,idim,im,iq) = os%zorb(is,idim,im)*os%phase(is, iq)
264 if(os%use_submesh)
then
265 do iorb = 1, os%norbs
266 call accel_write_buffer(os%buff_eorb(iq), ns, os%eorb_submesh(:, 1, iorb, iq), offset = (iorb - 1)*os%ldorbs_eorb)
269 do iorb = 1, os%norbs
271 os%eorb_mesh(:, iorb, 1, iq), offset = (iorb - 1)*os%ldorbs_eorb)
284 integer,
intent(in) :: dim
287 logical,
intent(in) :: spin_polarized
288 real(real64),
optional,
allocatable,
intent(in) :: vec_pot(:)
289 real(real64),
optional,
allocatable,
intent(in) :: vec_pot_var(:, :)
290 integer,
optional,
intent(in) :: kpt_max
292 integer :: iq, ikpoint
293 real(real64) :: kr, kpoint(dim), dx(dim)
294 integer :: inn, kpt_end
299 if(
present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
301 do iq = kpt%start, kpt_end
303 if(spin_polarized)
then
304 ikpoint = 1 + (iq - 1)/2
312 kpoint(1:dim) = kpoints%get_point(ikpoint)
314 if (os%nneighbors > 0)
then
315 do inn = 1, os%nneighbors
316 dx(1:dim) = os%V_ij(inn,1:dim)
317 kr = sum(kpoint(1:dim)*dx(1:dim))
318 if (
present(vec_pot))
then
319 if (
allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*dx(1:dim))
323 if (
present(vec_pot_var))
then
324 if (
allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, 1)*dx(1:dim))
328 os%phase_shift(inn, iq) =
exp(-
m_zi*kr)
345 type(
ions_t),
intent(in) :: ions
351 os%spec_index = os%spec%get_index()
352 do ja = 1, ions%natoms
353 if(ions%atom(ja)%species == os%spec)
then
354 os%spec_index = min(os%spec_index, ions%atom(ja)%species%get_index())
363#include "orbitalset_inc.F90"
366#include "complex.F90"
367#include "orbitalset_inc.F90"
double exp(double __x) __attribute__((__nothrow__
subroutine, public accel_release_buffer(this)
pure logical function, public accel_is_enabled()
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....
real(real64), parameter, public m_zero
complex(real64), parameter, public m_z0
complex(real64), parameter, public m_zi
real(real64), parameter, public m_one
integer pure function, public kpoints_number(this)
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.
subroutine, public zorbitalset_get_position_matrix_elem(os, ndim, psib, idir, dot)
subroutine, public orbitalset_update_phase_shift(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase shift for the intersite interaction.
subroutine orbitalset_set_species_index(os, ions)
Set the species index for a given orbital set.
subroutine, public orbitalset_init(this)
subroutine, public orbitalset_end(this)
subroutine, public dorbitalset_add_to_batch(os, ndim, psib, weight)
subroutine, public dorbitalset_get_position_matrix_elem(os, ndim, psib, idir, dot)
subroutine, public zorbitalset_add_to_batch(os, ndim, psib, weight)
subroutine, public orbitalset_update_phase(os, dim, kpt, kpoints, spin_polarized, vec_pot, vec_pot_var, kpt_max)
Build the phase correction to the global phase in case the orbital crosses the border of the simulato...
subroutine, public dorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public dorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public orbitalset_set_jln(this, jj, ll, nn)
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot, reduce)
subroutine, public zorbitalset_get_coefficients(os, ndim, psi, ik, has_phase, dot, reduce)
subroutine, public submesh_end(this)
Distribution of N instances over mpi_grpsize processes, for the local rank mpi_grprank....