Octopus
poisson_corrections.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
24 use global_oct_m
25 use iso_c_binding
28 use math_oct_m
30 use mesh_oct_m
34 use parser_oct_m
36 use space_oct_m
37
38 implicit none
39
40 private
41 public :: &
50
51 real(real64), parameter :: alpha_ = m_five
52
54 private
55 integer :: method
56 integer :: maxl
57 real(real64), allocatable :: phi(:, :)
58 real(real64), allocatable :: aux(:, :)
59 real(real64), allocatable :: gaussian(:)
60 end type poisson_corr_t
61
62 type(derivatives_t), pointer :: der_pointer
63 type(mesh_t), pointer :: mesh_pointer
64
65 integer, parameter :: &
66 CORR_MULTIPOLE = 1, &
67 corr_exact = 3
68
69contains
70
71 ! ---------------------------------------------------------
72 subroutine poisson_corrections_init(this, namespace, space, ml, mesh)
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
78
79 real(real64) :: alpha, beta, ylm, rr, xx(space%dim)
80 integer :: ip, ll, add_lm, lldfac, jj, mm
81
83
84 assert(space%dim == 3)
85
86 if (space%is_periodic()) then
87 call messages_not_implemented("Poisson boundary corrections for periodic systems", namespace=namespace)
88 end if
89
90 !%Variable PoissonSolverBoundaries
91 !%Type integer
92 !%Section Hamiltonian::Poisson
93 !%Default multipole
94 !%Description
95 !% For finite systems, some Poisson solvers (<tt>multigrid</tt>,
96 !% <tt>cg_corrected</tt>, and <tt>fft</tt> with <tt>PoissonFFTKernel = multipole_correction</tt>)
97 !% require the calculation of the
98 !% boundary conditions with an auxiliary method. This variable selects that method.
99 !%Option multipole 1
100 !% A multipole expansion of the density is used to approximate the potential on the boundaries.
101 !%Option exact 3
102 !% An exact integration of the Poisson equation is done over the boundaries. This option is
103 !% experimental, and not implemented for domain parallelization.
104 !%End
105 call parse_variable(namespace, 'PoissonSolverBoundaries', corr_multipole, this%method)
106
107 select case (this%method)
108 case (corr_multipole)
109 this%maxl = ml
110
111 add_lm = 0
112 do ll = 0, this%maxl
113 do mm = -ll, ll
114 add_lm = add_lm + 1
115 end do
116 end do
117
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))
121
122 alpha = alpha_ * mesh%spacing(1)
123 do ip = 1, mesh%np
124 call mesh_r(mesh, ip, rr, coords = xx)
125 this%gaussian(ip) = exp(-(rr/alpha)**2)
126 add_lm = 1
127 do ll = 0, this%maxl
128 lldfac = 1
129 do jj = 1, 2*ll+1, 2
130 lldfac = lldfac * jj
131 end do
132 beta = sqrt(m_pi)*2**(ll+3) / lldfac
133 do mm = -ll, ll
134 call ylmr_real(xx, ll, mm, ylm)
135 if (rr > m_epsilon) then
136 this%phi(ip, add_lm) = beta*isubl(ll, rr/alpha)*ylm/rr**(ll+1)
137 this%aux(ip, add_lm) = rr**ll*ylm
138 else
139 this%phi(ip, add_lm) = beta*ylm / alpha
140 if (ll == 0) then
141 this%aux(ip, add_lm) = ylm
142 else
143 this%aux(ip, add_lm) = m_zero
144 end if
145 end if
146 add_lm = add_lm + 1
147 end do
148 end do
149 end do
152 call messages_experimental('Exact Poisson solver boundaries', namespace=namespace)
153 if (mesh%parallel_in_domains) call messages_not_implemented('Exact Poisson solver boundaries with domain parallelization')
155 end select
156
159 contains
161 ! ---------------------------------------------------------
162 real(real64) function isubl( ll, xx)
163 integer, intent(in) :: ll
164 real(real64), intent(in) :: xx
165
166 ! no push_sub, called too frequently
167 isubl = m_half * gamma((ll + m_half)) * (m_one - loct_incomplete_gamma(ll+m_half, xx**2))
168
169 end function isubl
170
171 end subroutine poisson_corrections_init
172
173
174 ! ---------------------------------------------------------
175 subroutine poisson_corrections_end(this)
176 type(poisson_corr_t), intent(inout) :: this
177
179
180 select case (this%method)
181 case (corr_multipole)
182 safe_deallocate_a(this%phi)
183 safe_deallocate_a(this%aux)
184 safe_deallocate_a(this%gaussian)
185 case (corr_exact)
186 end select
187
189 end subroutine poisson_corrections_end
190
191
192 ! ---------------------------------------------------------
193 subroutine correct_rho(this, der, rho, rho_corrected, vh_correction)
194 type(poisson_corr_t), intent(in) :: this
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(:)
199
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(:)
204
205 push_sub(correct_rho)
206 call profiling_in("POISSON_CORRECT")
207
208 assert(ubound(vh_correction, dim = 1) == der%mesh%np_part)
209
210 select case (this%method)
211 case (corr_multipole)
212
213 safe_allocate(mult(1:(this%maxl+1)**2))
214 call get_multipoles(this, der%mesh, rho, this%maxl, mult)
215
216 alpha = alpha_ * der%mesh%spacing(1)
217
218 safe_allocate(betal(1:(this%maxl+1)**2))
219 add_lm = 1
220 do ll = 0, this%maxl
221 do mm = -ll, ll
222 lldfac = 1
223 do jj = 1, 2*ll+1, 2
224 lldfac = lldfac*jj
225 end do
226 betal(add_lm) = (2**(ll + 2)) / (alpha**(2*ll+3) * sqrt(m_pi) * lldfac)
227 add_lm = add_lm + 1
228 end do
229 end do
230
231 call lalg_copy(der%mesh%np, rho, rho_corrected)
232 vh_correction = m_zero
233 add_lm = 1
234 do ll = 0, this%maxl
235 do mm = -ll, ll
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)
238 end do
239 call lalg_axpy(der%mesh%np, mult(add_lm), this%phi(:, add_lm), vh_correction)
240 add_lm = add_lm + 1
241 end do
242 end do
243
244 safe_deallocate_a(mult)
245 safe_deallocate_a(betal)
246
247 case (corr_exact)
248
249 do ip = 1, der%mesh%np
250 vh_correction(ip) = m_zero
251 end do
252
253 do ip = der%mesh%np + 1, der%mesh%np_part
254 vv = m_zero
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
258 end do
259 vh_correction(ip) = der%mesh%volume_element*vv
260 end do
261
262 call dderivatives_lapl(der, vh_correction, rho_corrected, set_bc = .false.)
263
264 do ip = 1, der%mesh%np
265 rho_corrected(ip) = rho(ip) + m_one/(m_four*m_pi)*rho_corrected(ip)
266 end do
267
268 end select
269
270 call profiling_out("POISSON_CORRECT")
271 pop_sub(correct_rho)
272 end subroutine correct_rho
273
274
275 ! ---------------------------------------------------------
276 subroutine get_multipoles(this, mesh, rho, ml, mult)
277 type(poisson_corr_t), intent(in) :: this
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)
282
283 integer :: add_lm, ll, mm
284
285 push_sub(get_multipoles)
286
287 mult(:) = m_zero
288 add_lm = 1
289 do ll = 0, ml
290 do mm = -ll, ll
291 mult(add_lm) = dmf_dotp(mesh, rho, this%aux(:, add_lm))
292 add_lm = add_lm + 1
293 end do
294 end do
295
296 pop_sub(get_multipoles)
297 end subroutine get_multipoles
298
299
300 ! ---------------------------------------------------------
301 real(real64) function internal_dotp(xx, yy) result(res)
302 real(real64), intent(in) :: xx(:)
303 real(real64), intent(in) :: yy(:)
304
305 push_sub(internal_dotp)
306
307 res = dmf_dotp(mesh_pointer, xx, yy)
308 pop_sub(internal_dotp)
309 end function internal_dotp
310
311
312 ! ---------------------------------------------------------
313 subroutine poisson_boundary_conditions(this, mesh, rho, pot)
314 type(poisson_corr_t), intent(in) :: this
315 type(mesh_t), intent(in) :: mesh
316 real(real64), intent(in) :: rho(:)
317 real(real64), intent(inout) :: pot(:)
318
319 integer :: ip, add_lm, ll, mm, bp_lower
320 real(real64) :: xx(mesh%box%dim), rr, s1, sa
321 real(real64), allocatable :: mult(:)
322
324
325 assert(mesh%box%dim == 3)
326
327 safe_allocate(mult(1:(this%maxl+1)**2))
328
329 call get_multipoles(this, mesh, rho, this%maxl, mult)
330
331 bp_lower = mesh%np + 1
332 if (mesh%parallel_in_domains) then
333 bp_lower = bp_lower + mesh%pv%np_ghost
334 end if
335
336 pot(bp_lower:mesh%np_part) = m_zero
337 do ip = bp_lower, mesh%np_part ! boundary conditions
338 call mesh_r(mesh, ip, rr, coords = xx)
339 add_lm = 1
340 do ll = 0, this%maxl
341 s1 = m_four*m_pi/((m_two*ll + m_one)*rr**(ll + 1))
342 do mm = -ll, ll
343 call ylmr_real(xx, ll, mm, sa)
344 pot(ip) = pot(ip) + sa * mult(add_lm) * s1
345 add_lm = add_lm+1
346 end do
347 end do
348 end do
349
350 safe_deallocate_a(mult)
352 end subroutine poisson_boundary_conditions
353
355
356!! Local Variables:
357!! mode: f90
358!! coding: utf-8
359!! End:
Define which routines can be seen from the outside.
Definition: loct_math.F90:148
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
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_five
Definition: global.F90:196
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public ylmr_real(xx, li, mi, ylm)
This is a Numerical Recipes-based subroutine computes real spherical harmonics ylm at position (x,...
Definition: math.F90:378
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:342
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
This module defines non-local operators.
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)
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)