Octopus
poisson_psolver.F90
Go to the documentation of this file.
1!! Copyright (C) 2013 J. Alberdi-Rodriguez
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
23 use cube_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
30 use mesh_oct_m
32 use math_oct_m
34 use parser_oct_m
36 use space_oct_m
38 ! Support for PSolver from BigDFT
39#ifdef HAVE_PSOLVER
41 use poisson_solver
42 use dictionaries, dict_set => set
43 use yaml_output, only: yaml_map
44 use wrapper_mpi, only: mpi_environment, mpi_environment_set
45 use at_domain
46#endif
47
48
49 implicit none
50
51 private
52 public :: &
60
62 private
63 type(fourier_space_op_t) :: coulb
64#ifdef HAVE_PSOLVER
65 type(coulomb_operator) :: kernel
66 type(dictionary), pointer :: inputs
67#endif
68
80 character(len = 1) :: geocode = "F"
90 character(len = 1), public :: datacode = "G"
91
92 integer :: isf_order
93 integer :: localnscatterarr(5)
94 real(real64) :: offset
95 end type poisson_psolver_t
96
97#ifdef HAVE_PSOLVER
98 logical, save :: flib_initialized = .false.
99#endif
100
101 real(real64), parameter :: TOL_VANISHING_Q = 1e-6_real64
102
103contains
104
105 !-----------------------------------------------------------------
106 subroutine poisson_psolver_init(this, namespace, space, cube, mu, qq, force_isolated)
107 type(poisson_psolver_t), intent(out) :: this
108 type(namespace_t), intent(in) :: namespace
109 class(space_t), intent(in) :: space
110 type(cube_t), intent(inout) :: cube
111 real(real64), intent(in) :: mu
112 real(real64), intent(in) :: qq(:)
113 logical, optional, intent(in) :: force_isolated
114
115 logical data_is_parallel
116#ifdef HAVE_PSOLVER
117 real(real64) :: alpha, beta, gamma
118 real(real64) :: modq2
119 type(mpi_environment) mpi_env
120 type(domain) :: dom
121#endif
122
123 push_sub(poisson_psolver_init)
124
125#ifdef HAVE_PSOLVER
126 if (.not. flib_initialized) then
127 call f_lib_initialize()
128 flib_initialized = .true.
129 end if
130 call dict_init(this%inputs)
131#endif
132
133 if (optional_default(force_isolated, .false.)) then
134 this%geocode = "F"
135 else
136 select case (space%periodic_dim)
137 case (0)
138 ! Free BC
139 this%geocode = "F"
140 case (1)
141 ! Wire BC
142 this%geocode = "W"
143 call messages_not_implemented("PSolver support for 1D periodic boundary conditions", namespace=namespace)
144 case (2)
145 ! Surface BC
146 this%geocode = "S"
147 call messages_not_implemented("PSolver support for 2D periodic boundary conditions", namespace=namespace)
148 case (3)
149 ! Periodic BC
150 this%geocode = "P"
151 call messages_experimental("PSolver support for 3D periodic boundary conditions", namespace=namespace)
152 end select
153 end if
154
155#ifdef HAVE_PSOLVER
156 !Verbosity switch
157 call dict_set(this%inputs//'setup'//'verbose', .false.)
158 !Order of the Interpolating Scaling Function family
159 call dict_set(this%inputs//'kernel'//'isf_order', 16)
160 !Mu screening parameter
161 call dict_set(this%inputs//'kernel'//'screening', mu)
162 !Calculation of the stress tensor
163 call dict_set(this%inputs//'kernel'//'stress_tensor', .false.)
164#else
165 this%isf_order = 16
166#endif
167
168 !%Variable PoissonSolverPSolverParallelData
169 !%Type logical
170 !%Section Hamiltonian::Poisson::PSolver
171 !%Default yes
172 !%Description
173 !% Indicates whether data is partitioned within the PSolver library.
174 !% If data is distributed among processes, Octopus uses parallel data-structures
175 !% and, thus, less memory.
176 !% If "yes", data is parallelized. The <i>z</i>-axis of the input vector
177 !% is split among the MPI processes.
178 !% If "no", entire input and output vector is saved in all the MPI processes.
179 !% If k-points parallelization is used, "no" must be selected.
180 !%End
181 call parse_variable(namespace, 'PoissonSolverPSolverParallelData', .true., data_is_parallel)
182 if (.not. cube%parallel_in_domains) then
183 data_is_parallel = .false.
184 end if
185
186 call messages_obsolete_variable(namespace, 'PoissonSolverISFParallelData', 'PoissonSolverPSolverParallelData')
188#ifdef HAVE_PSOLVER
189 if (data_is_parallel) then
190 call dict_set(this%inputs//'setup'//'global_data', .false.)
191 else
192 call dict_set(this%inputs//'setup'//'global_data', .true.)
193 end if
194#else
195 if (data_is_parallel) then
196 this%datacode = "D"
197 else
198 this%datacode = "G"
199 end if
200#endif
201
202#ifdef HAVE_PSOLVER
203 call dict_set(this%inputs//'setup'//'verbose', debug%info)
204
205 alpha = cube%latt%alpha*m_pi/(180.0_real64)
206 beta = cube%latt%beta*m_pi/(180.0_real64)
207 gamma = cube%latt%gamma*m_pi/(180.0_real64)
208
209 ! Previously, pkernel_init set the communicator used within PSolver to comm_world.
210 ! This can be overwritten by passing an optional argument of type(mpi_environment)
211 ! to pkernel_init(). This data type is defined within the wrapper_MPI module of
212 ! the Futile library. Futile is a prerequisit for PSolver.
213
214 ! TODO: check that cube%mpi_grp corresponds to the correct mpi group, when parallelizing, e.g., over systems !!
215
216 call mpi_environment_set(mpi_env, cube%mpi_grp%rank, cube%mpi_grp%size, cube%mpi_grp%comm%MPI_VAL, cube%mpi_grp%size)
217
218 dom=domain_new(units=atomic_units, bc=geocode_to_bc_enum(this%geocode),&
219 alpha_bc=alpha, beta_ac=beta, gamma_ab=gamma, acell=cube%rs_n_global*cube%spacing)
220
221 this%kernel = pkernel_init(cube%mpi_grp%rank, cube%mpi_grp%size, this%inputs, dom, cube%rs_n_global, &
222 cube%spacing, alpha_bc = alpha, beta_ac = beta, gamma_ab = gamma, mpi_env = mpi_env)
223
224 call pkernel_set(this%kernel, verbose=debug%info)
225
226 !G=0 component
227 modq2 = sum(qq(1:space%periodic_dim)**2)
228 if (modq2 > m_epsilon) then
229 this%offset = m_one/modq2
230 else
231 this%offset = m_zero
232 end if
233
234 !Screened coulomb potential (erfc function)
235 if (mu > m_epsilon) then
236 if (modq2 > m_epsilon) then
237 this%offset = this%offset*(-expm1(-modq2/((m_two*mu)**2)))
238 else
239 !Analytical limit of 1/|q|^2*(1-exp(-|q|^2/4mu^2))
240 this%offset = m_one/((m_two*mu)**2)
241 end if
242 end if
243 this%offset = this%offset*m_four*m_pi
244#endif
245
246 pop_sub(poisson_psolver_init)
247 end subroutine poisson_psolver_init
248
249 !-----------------------------------------------------------------
250 subroutine poisson_psolver_end(this)
251 type(poisson_psolver_t), intent(inout) :: this
252
253 push_sub(poisson_psolver_end)
254
255#ifdef HAVE_PSOLVER
256 call pkernel_free(this%kernel)
257 call f_lib_finalize()
258 call dict_free(this%inputs)
259#endif
260
261 pop_sub(poisson_psolver_end)
262 end subroutine poisson_psolver_end
263
264 !-----------------------------------------------------------------
265 subroutine poisson_psolver_reinit(this, space, cube, coulb, qq_in)
266 type(poisson_psolver_t), intent(inout) :: this
267 class(space_t), intent(in) :: space
268 type(cube_t), intent(inout) :: cube
269 type(fourier_space_op_t), intent(inout) :: coulb
270 real(real64), intent(in) :: qq_in(:)
271
272 real(real64) :: alpha, beta, gamma
273 real(real64) :: qq_abs(1:space%dim), qq(1:space%dim)
274 real(real64) :: modq2
275 integer :: idim
276#ifdef HAVE_PSOLVER
277 type(domain) :: dom
278#endif
279 push_sub(poisson_psolver_reinit)
280
281 ! We only support short range only or long-range only at the moment
282 if(coulb%mu > m_epsilon) then
283 assert(coulb%alpha < m_epsilon .or. coulb%beta < m_epsilon)
284 end if
285
286#ifdef HAVE_PSOLVER
287 call pkernel_free(this%kernel)
288#endif
289
290 !We might change the cell angles
291 alpha = cube%latt%alpha*m_pi/(180.0_real64)
292 beta = cube%latt%beta*m_pi/(180.0_real64)
293 gamma = cube%latt%gamma*m_pi/(180.0_real64)
294
295#ifdef HAVE_PSOLVER
296 call dict_set(this%inputs//'kernel'//'screening', coulb%mu)
297
298 dom=domain_new(units=atomic_units, bc=geocode_to_bc_enum(this%geocode),&
299 alpha_bc=alpha, beta_ac=beta, gamma_ab=gamma, acell=cube%rs_n_global*cube%spacing)
300
301 this%kernel = pkernel_init(cube%mpi_grp%rank, cube%mpi_grp%size, this%inputs,&
302 dom,cube%rs_n_global, cube%spacing, &
303 alpha_bc = alpha, beta_ac = beta, gamma_ab = gamma)
304
305 call pkernel_set(this%kernel, verbose=debug%info)
306#endif
307
308 !G=0 component
309 !We remove potential umklapp
310 do idim = 1, space%periodic_dim
311 qq(idim) = qq_in(idim) - anint(qq_in(idim) + m_half*1e-8_real64)
312 end do
313 qq(space%periodic_dim + 1:space%dim) = m_zero
314 call kpoints_to_absolute(cube%latt, qq, qq_abs)
315 modq2 = norm2(qq_abs)
316 if (modq2 > tol_vanishing_q) then
317 this%offset = m_four*m_pi/modq2
318 else
319 this%offset = coulb%singularity
320 end if
321
322 !Screened coulomb potential (erfc function)
323 if (coulb%mu > m_epsilon) then
324 if (modq2 > tol_vanishing_q) then
325 this%offset = this%offset*(-expm1(-modq2/((m_two*coulb%mu)**2)))
326 else
327 !Analytical limit of 1/|q|^2*(1-exp(-|q|^2/4mu^2))
328 this%offset = m_four*m_pi/((m_two*coulb%mu)**2)
329 end if
330 end if
331
333 end subroutine poisson_psolver_reinit
334
335 !-----------------------------------------------------------------
336 subroutine poisson_psolver_parallel_solve(this, mesh, cube, pot, rho, mesh_cube_map)
337 type(poisson_psolver_t), intent(in), target :: this
338 type(mesh_t), intent(in) :: mesh
339 type(cube_t), intent(in) :: cube
340 real(real64), intent(out) :: pot(:)
341 real(real64), intent(in) :: rho(:)
342 type(mesh_cube_parallel_map_t), intent(in) :: mesh_cube_map
343
344 type(cube_function_t) :: cf
345
346#ifdef HAVE_PSOLVER
347 type(coulomb_operator), pointer :: kernel_pointer
348 real(real64) :: hartree_energy
349 real(real64) :: offset
350 ! To be used only in the periodic case, ignored for other boundary conditions.
351#endif
352
356 real(real64), allocatable :: pot_ion(:,:,:)
357 character(len=3) :: quiet
358
360
361 call dcube_function_alloc_rs(cube, cf)
362 call dmesh_to_cube_parallel(mesh, rho, cube, cf, mesh_cube_map)
363 safe_allocate(pot_ion(1:cube%rs_n(1),1:cube%rs_n(2),1:cube%rs_n(3)))
364
365 quiet = "YES"
366 if (debug%info) quiet = "NO"
367
368#ifdef HAVE_PSOLVER
369 !The offset is the integral over space of the potential
370 !this%offset is the G=0 component of the (screened) Coulomb potential
371 !The G=0 component of the Hartree therefore needs to be
372 ! multiplied by the G=0 component of the density
373 if (this%offset > m_epsilon) then
374 offset = this%offset*dmf_integrate(mesh,rho)
375 end if
376#endif
377
378 call profiling_in("PSOLVER_LIBRARY")
379#ifdef HAVE_PSOLVER
380 kernel_pointer => this%kernel
381 call h_potential(this%datacode, kernel_pointer, cf%dRS, pot_ion, hartree_energy, offset, .false., quiet = quiet)
382#endif
383 call profiling_out("PSOLVER_LIBRARY")
384 safe_deallocate_a(pot_ion)
385
386 call dcube_to_mesh_parallel(cube, cf, mesh, pot, mesh_cube_map)
388 call dcube_function_free_rs(cube, cf)
389
391 end subroutine poisson_psolver_parallel_solve
392
393 !-----------------------------------------------------------------
394 subroutine poisson_psolver_global_solve(this, mesh, cube, pot, rho, sm)
395 type(poisson_psolver_t), intent(in), target :: this
396 type(mesh_t), intent(in) :: mesh
397 type(cube_t), intent(in) :: cube
398 real(real64), contiguous, intent(out) :: pot(:)
399 real(real64), contiguous, intent(in) :: rho(:)
400 type(submesh_t), optional, intent(in) :: sm
401
402 character(len=3) :: quiet
403 type(cube_function_t) :: cf
404
405#ifdef HAVE_PSOLVER
406 type(coulomb_operator), pointer :: kernel_pointer
407 real(real64) :: hartree_energy
408 real(real64) :: offset
409 ! To be used only in the periodic case, ignored for other boundary conditions.
410#endif
411
416 real(real64), allocatable :: pot_ion(:,:,:)
417
419
420 call dcube_function_alloc_rs(cube, cf)
421
422 if (present(sm)) then
423 call dsubmesh_to_cube(sm, rho, cube, cf)
424 else
425 call dmesh_to_cube(mesh, rho, cube, cf)
426 end if
427
428 safe_allocate(pot_ion(1:cube%rs_n(1), 1:cube%rs_n(2), 1:cube%rs_n(3)))
429
430 quiet = "YES"
431 if (debug%info) quiet = "NO"
432
433#ifdef HAVE_PSOLVER
434 !The offset is the integral over space of the potential
435 !this%offset is the G=0 component of the (screened) Coulomb potential
436 !The G=0 component of the Hartree therefore needs to be
437 ! multiplied by the G=0 component of the density
438 if (this%offset > m_epsilon) then
439 offset = this%offset*dmf_integrate(mesh, rho)
440 end if
441#endif
442
443 call profiling_in("PSOLVER_LIBRARY")
444
445#ifdef HAVE_PSOLVER
446 kernel_pointer => this%kernel
447 call h_potential(this%datacode, kernel_pointer, &
448 cf%dRS, pot_ion, hartree_energy, offset, .false., &
449 quiet = quiet) !optional argument
450#endif
451 safe_deallocate_a(pot_ion)
452
453 if (present(sm)) then
454 call dcube_to_submesh(cube, cf, sm, pot)
455 else
456 call dcube_to_mesh(cube, cf, mesh, pot)
457 end if
458
459 call dcube_function_free_rs(cube, cf)
460
462 end subroutine poisson_psolver_global_solve
463
464 subroutine poisson_psolver_get_dims(this, cube)
465 type(poisson_psolver_t), intent(inout) :: this
466 type(cube_t), intent(inout) :: cube
467
468#ifdef HAVE_PSOLVER
469
477 integer :: n3d
480 integer :: n3p
483 integer :: n3pi
489 integer :: i3xcsh
497 integer :: i3s
498
500 logical :: use_gradient
502 logical :: use_wb_corr
503#endif
504
506
507 ! Get the dimensions of the cube
508#ifdef HAVE_PSOLVER
509 use_gradient = .false.
510 use_wb_corr = .false.
511 call ps_dim4allocation(this%geocode, this%datacode, cube%mpi_grp%rank, cube%mpi_grp%size, &
512 cube%rs_n_global(1), cube%rs_n_global(2), cube%rs_n_global(3), &
513 use_gradient, use_wb_corr, &
514 0, n3d, n3p, n3pi, i3xcsh, i3s)
515 this%localnscatterarr(:) = (/ n3d, n3p, n3pi, i3xcsh, i3s /)
516#endif
517 cube%rs_n(1:2) = cube%rs_n_global(1:2)
518 cube%rs_n(3) = this%localnscatterarr(1)
519 cube%rs_istart(1:2) = 1
520 cube%rs_istart(3) = this%localnscatterarr(5)
521
522 !! With PSolver we don`t care about the Fourier space and its dimensions
523 !! We`ll put as in RS
524 cube%fs_n_global(1) = cube%rs_n_global(1)
525 cube%fs_n_global(2) = cube%rs_n_global(2)
526 cube%fs_n_global(3) = cube%rs_n_global(3)
527 cube%fs_n(1:2) = cube%rs_n_global(1:2)
528 cube%fs_n(3) = this%localnscatterarr(1)
529 cube%fs_istart(1:2) = 1
530 cube%fs_istart(3) = this%localnscatterarr(5)
531
533 end subroutine poisson_psolver_get_dims
534
535end module poisson_psolver_oct_m
536
537!! Local Variables:
538!! mode: f90
539!! coding: utf-8
540!! End:
subroutine, public dmesh_to_cube(mesh, mf, cube, cf)
Convert a function from the mesh to the cube.
subroutine, public dcube_to_submesh(cube, cf, sm, mf)
subroutine, public dcube_to_mesh(cube, cf, mesh, mf)
Convert a function from the cube to the mesh.
subroutine, public dcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
subroutine, public dcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public dsubmesh_to_cube(sm, mf, cube, cf)
The next two subroutines convert a function between a submesh and the cube.
subroutine, public dmesh_to_cube_parallel(mesh, mf, cube, cf, map)
The next two subroutines convert a function between the normal mesh and the cube in parallel.
subroutine, public dcube_to_mesh_parallel(cube, cf, mesh, mf, map)
type(debug_t), save, public debug
Definition: debug.F90:158
real(real64), parameter, public m_two
Definition: global.F90:202
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public m_four
Definition: global.F90:204
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:198
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1137
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public poisson_psolver_global_solve(this, mesh, cube, pot, rho, sm)
subroutine, public poisson_psolver_parallel_solve(this, mesh, cube, pot, rho, mesh_cube_map)
subroutine, public poisson_psolver_end(this)
subroutine, public poisson_psolver_reinit(this, space, cube, coulb, qq_in)
subroutine, public poisson_psolver_get_dims(this, cube)
subroutine, public poisson_psolver_init(this, namespace, space, cube, mu, qq, force_isolated)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
Describes mesh distribution to nodes.
Definition: mesh.F90:187
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:175
int true(void)