73 type(submesh_t) :: sphere
74 complex(real64),
allocatable :: phase(:,:)
77 real(real64) :: Ubar, Jbar
81 real(real64),
allocatable :: V_ij(:,:)
82 real(real64),
allocatable :: coulomb_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(:,:,:,:)
95 type(accel_mem_t) :: dbuff_orb, zbuff_orb
96 type(accel_mem_t),
allocatable :: buff_eorb (:)
98 logical :: use_submesh
100 type(poisson_t) :: poisson
106 type(orbitalset_t),
intent(inout) :: this
132 type(orbitalset_t),
intent(inout) :: this
138 safe_deallocate_a(this%phase)
139 safe_deallocate_a(this%dorb)
140 safe_deallocate_a(this%zorb)
141 safe_deallocate_a(this%eorb_submesh)
142 safe_deallocate_a(this%eorb_mesh)
146 safe_deallocate_a(this%V_ij)
147 safe_deallocate_a(this%coulomb_IIJJ)
148 safe_deallocate_a(this%map_os)
149 safe_deallocate_a(this%phase_shift)
154 if(
allocated(this%buff_eorb))
then
155 do ik = lbound(this%buff_eorb, dim=1), ubound(this%buff_eorb, dim=1)
158 safe_deallocate_a(this%buff_eorb)
167 real(real64),
intent(in) :: jj
168 integer,
intent(in) :: ll, nn
183 integer,
intent(in) :: dim
186 logical,
intent(in) :: spin_polarized
187 real(real64),
optional,
allocatable,
intent(in) :: vec_pot(:)
188 real(real64),
optional,
allocatable,
intent(in) :: vec_pot_var(:, :)
189 integer,
optional,
intent(in) :: kpt_max
191 integer :: ns, iq, is, ikpoint, im, idim, kpt_end
192 real(real64) :: kr, kpoint(1:dim), dx(1:dim)
200 if (
present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
202 do iq = kpt%start, kpt_end
204 if (spin_polarized)
then
205 ikpoint = 1 + (iq - 1)/2
213 kpoint(1:dim) = kpoints%get_point(ikpoint)
219 dx = os%sphere%rel_x(1:dim, is) - os%sphere%mesh%x(os%sphere%map(is), 1:dim) + os%sphere%center(1:dim)
220 kr = dot_product(kpoint, dx)
221 if (
present(vec_pot))
then
222 if (
allocated(vec_pot)) kr = kr + dot_product(vec_pot, dx)
225 if (
present(vec_pot_var))
then
226 if (
allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, os%sphere%map(is)) &
227 *(os%sphere%rel_x(1:dim, is)+os%sphere%center))
230 os%phase(is, iq) =
exp(
m_zi*kr)
233 if (.not. os%use_submesh)
then
237 os%eorb_mesh(:, im, idim, iq) =
m_z0
239 os%eorb_mesh(os%sphere%map(is),im,idim,iq) = os%eorb_mesh(os%sphere%map(is),im,idim,iq) &
240 + os%zorb(is, idim, im) * os%phase(is, iq)
250 os%eorb_submesh(is,idim,im,iq) = os%zorb(is,idim,im)*os%phase(is, iq)
259 if(os%use_submesh)
then
260 do iorb = 1, os%norbs
261 call accel_write_buffer(os%buff_eorb(iq), ns, os%eorb_submesh(:, 1, iorb, iq), offset = (iorb - 1)*os%ldorbs)
264 do iorb = 1, os%norbs
266 os%eorb_mesh(:, iorb, 1, iq), offset = (iorb - 1)*os%ldorbs)
279 integer,
intent(in) :: dim
282 logical,
intent(in) :: spin_polarized
283 real(real64),
optional,
allocatable,
intent(in) :: vec_pot(:)
284 real(real64),
optional,
allocatable,
intent(in) :: vec_pot_var(:, :)
285 integer,
optional,
intent(in) :: kpt_max
287 integer :: iq, ikpoint
288 real(real64) :: kr, kpoint(dim), dx(dim)
289 integer :: inn, kpt_end
294 if(
present(kpt_max)) kpt_end = min(kpt_max, kpt_end)
296 do iq = kpt%start, kpt_end
298 if(spin_polarized)
then
299 ikpoint = 1 + (iq - 1)/2
307 kpoint(1:dim) = kpoints%get_point(ikpoint)
309 if (os%nneighbors > 0)
then
310 do inn = 1, os%nneighbors
311 dx(1:dim) = os%V_ij(inn,1:dim)
312 kr = sum(kpoint(1:dim)*dx(1:dim))
313 if (
present(vec_pot))
then
314 if (
allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*dx(1:dim))
318 if (
present(vec_pot_var))
then
319 if (
allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, 1)*dx(1:dim))
323 os%phase_shift(inn, iq) =
exp(-
m_zi*kr)
334#include "orbitalset_inc.F90"
337#include "complex.F90"
338#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, public dorbitalset_get_coeff_batch(os, ndim, psib, dot)
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_accel(os, ndim, psib, dot)
subroutine, public orbitalset_set_jln(this, jj, ll, nn)
subroutine, public zorbitalset_get_coeff_batch_accel(os, ndim, psib, dot)
subroutine, public zorbitalset_get_coeff_batch(os, ndim, psib, dot)
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....