51 real(real64),
parameter :: alpha_ =
m_five
57 real(real64),
allocatable :: phi(:, :)
58 real(real64),
allocatable :: aux(:, :)
59 real(real64),
allocatable :: gaussian(:)
62 type(derivatives_t),
pointer :: der_pointer
63 type(mesh_t),
pointer :: mesh_pointer
65 integer,
parameter :: &
73 type(poisson_corr_t),
intent(out) :: this
74 type(namespace_t),
intent(in) :: namespace
75 class(space_t),
intent(in) :: space
76 integer,
intent(in) :: ml
77 type(mesh_t),
intent(in) :: mesh
79 real(real64) :: alpha, beta, ylm, rr, xx(space%dim)
80 integer :: ip, ll, add_lm, lldfac, jj, mm
84 assert(space%dim == 3)
86 if (space%is_periodic())
then
105 call parse_variable(namespace,
'PoissonSolverBoundaries', corr_multipole, this%method)
107 select case (this%method)
108 case (corr_multipole)
118 safe_allocate(this%phi(1:mesh%np, 1:add_lm))
119 safe_allocate(this%aux(1:mesh%np, 1:add_lm))
120 safe_allocate(this%gaussian(1:mesh%np))
122 alpha = alpha_ * mesh%spacing(1)
124 call mesh_r(mesh, ip, rr, coords = xx)
125 this%gaussian(ip) =
exp(-(rr/alpha)**2)
132 beta =
sqrt(
m_pi)*2**(ll+3) / lldfac
136 this%phi(ip, add_lm) = beta*
isubl(ll, rr/alpha)*ylm/rr**(ll+1)
137 this%aux(ip, add_lm) = rr**ll*ylm
139 this%phi(ip, add_lm) = beta*ylm / alpha
141 this%aux(ip, add_lm) = ylm
143 this%aux(ip, add_lm) =
m_zero
162 real(real64) function
isubl( ll, xx)
163 integer,
intent(in) :: ll
164 real(real64),
intent(in) :: xx
180 select case (this%method)
182 safe_deallocate_a(this%phi)
183 safe_deallocate_a(this%aux)
184 safe_deallocate_a(this%gaussian)
193 subroutine correct_rho(this, der, rho, rho_corrected, vh_correction)
195 type(derivatives_t),
intent(in) :: der
196 real(real64),
contiguous,
intent(in) :: rho(:)
197 real(real64),
contiguous,
intent(out) :: rho_corrected(:)
198 real(real64),
contiguous,
intent(out) :: vh_correction(:)
200 integer :: ip, add_lm, ll, mm, lldfac, jj, ip2
201 real(real64) :: alpha, vv, rr
202 real(real64),
allocatable :: mult(:)
203 real(real64),
allocatable :: betal(:)
206 call profiling_in(
"POISSON_CORRECT")
208 assert(ubound(vh_correction, dim = 1) == der%mesh%np_part)
210 select case (this%method)
213 safe_allocate(mult(1:(this%maxl+1)**2))
216 alpha =
alpha_ * der%mesh%spacing(1)
218 safe_allocate(betal(1:(this%maxl+1)**2))
226 betal(add_lm) = (2**(ll + 2)) / (alpha**(2*ll+3) *
sqrt(m_pi) * lldfac)
231 call lalg_copy(der%mesh%np, rho, rho_corrected)
232 vh_correction = m_zero
236 do ip = 1, der%mesh%np
237 rho_corrected(ip) = rho_corrected(ip) - mult(add_lm)*betal(add_lm)*this%aux(ip, add_lm)*this%gaussian(ip)
239 call lalg_axpy(der%mesh%np, mult(add_lm), this%phi(:, add_lm), vh_correction)
244 safe_deallocate_a(mult)
245 safe_deallocate_a(betal)
249 do ip = 1, der%mesh%np
250 vh_correction(ip) = m_zero
253 do ip = der%mesh%np + 1, der%mesh%np_part
255 do ip2 = 1, der%mesh%np
256 rr = norm2((der%mesh%x(1:der%dim, ip) - der%mesh%x(1:der%dim, ip2)))
257 vv = vv + rho(ip2)/rr
259 vh_correction(ip) = der%mesh%volume_element*vv
262 call dderivatives_lapl(der, vh_correction, rho_corrected, set_bc = .false.)
264 do ip = 1, der%mesh%np
265 rho_corrected(ip) = rho(ip) + m_one/(m_four*m_pi)*rho_corrected(ip)
270 call profiling_out(
"POISSON_CORRECT")
278 type(mesh_t),
intent(in) :: mesh
279 real(real64),
intent(in) :: rho(:)
280 integer,
intent(in) :: ml
281 real(real64),
intent(out) :: mult((ml+1)**2)
283 integer :: add_lm, ll, mm
291 mult(add_lm) = dmf_dotp(mesh, rho, this%aux(:, add_lm))
302 real(real64),
intent(in) :: xx(:)
303 real(real64),
intent(in) :: yy(:)
315 type(mesh_t),
intent(in) :: mesh
316 real(real64),
intent(in) :: rho(:)
317 real(real64),
intent(inout) :: pot(:)
319 integer :: ip, add_lm, ll, mm, bp_lower
320 real(real64) :: xx(mesh%box%dim), rr, s1, sa
321 real(real64),
allocatable :: mult(:)
325 assert(mesh%box%dim == 3)
327 safe_allocate(mult(1:(this%maxl+1)**2))
331 bp_lower = mesh%np + 1
332 if (mesh%parallel_in_domains)
then
333 bp_lower = bp_lower + mesh%pv%np_ghost
336 pot(bp_lower:mesh%np_part) = m_zero
337 do ip = bp_lower, mesh%np_part
338 call mesh_r(mesh, ip, rr, coords = xx)
341 s1 = m_four*m_pi/((m_two*ll + m_one)*rr**(ll + 1))
343 call ylmr_real(xx, ll, mm, sa)
344 pot(ip) = pot(ip) + sa * mult(add_lm) * s1
350 safe_deallocate_a(mult)
Define which routines can be seen from the outside.
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_five
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public ylmr_real(xx, li, mi, ylm)
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x,...
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
subroutine, public messages_not_implemented(feature, namespace)
subroutine, public messages_experimental(name, namespace)
This module defines non-local operators.
integer, parameter corr_multipole
subroutine get_multipoles(this, mesh, rho, ml, mult)
subroutine, public poisson_corrections_end(this)
real(real64) function, public internal_dotp(xx, yy)
subroutine, public poisson_boundary_conditions(this, mesh, rho, pot)
subroutine, public poisson_corrections_init(this, namespace, space, ml, mesh)
integer, parameter corr_exact
real(real64), parameter alpha_
type(mesh_t), pointer, public mesh_pointer
subroutine, public correct_rho(this, der, rho, rho_corrected, vh_correction)
type(derivatives_t), pointer, public der_pointer
real(real64) function isubl(ll, xx)