27 use,
intrinsic :: iso_fortran_env
42 use dictionaries, dict_set => set
43 use yaml_output,
only: yaml_map
44 use wrapper_mpi,
only: mpi_environment, mpi_environment_set
63 type(fourier_space_op_t) :: coulb
65 type(coulomb_operator) :: kernel
66 type(dictionary),
pointer :: inputs
80 character(len = 1) :: geocode =
"F"
90 character(len = 1),
public :: datacode =
"G"
93 integer :: localnscatterarr(5)
94 real(real64) :: offset
98 logical,
save :: flib_initialized = .false.
101 real(real64),
parameter :: TOL_VANISHING_Q = 1e-6_real64
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
115 logical data_is_parallel
117 real(real64) :: alpha, beta, gamma
118 real(real64) :: modq2
119 type(mpi_environment) mpi_env
126 if (.not. flib_initialized)
then
127 call f_lib_initialize()
128 flib_initialized = .
true.
130 call dict_init(this%inputs)
136 select case (space%periodic_dim)
151 call messages_experimental(
"PSolver support for 3D periodic boundary conditions", namespace=namespace)
157 call dict_set(this%inputs//
'setup'//
'verbose', .false.)
159 call dict_set(this%inputs//
'kernel'//
'isf_order', 16)
161 call dict_set(this%inputs//
'kernel'//
'screening', mu)
163 call dict_set(this%inputs//
'kernel'//
'stress_tensor', .false.)
181 call parse_variable(namespace,
'PoissonSolverPSolverParallelData', .
true., data_is_parallel)
182 if (.not. cube%parallel_in_domains)
then
183 data_is_parallel = .false.
189 if (data_is_parallel)
then
190 call dict_set(this%inputs//
'setup'//
'global_data', .false.)
192 call dict_set(this%inputs//
'setup'//
'global_data', .
true.)
195 if (data_is_parallel)
then
203 call dict_set(this%inputs//
'setup'//
'verbose',
debug%info)
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)
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)
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)
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)
224 call pkernel_set(this%kernel, verbose=
debug%info)
227 modq2 = sum(qq(1:space%periodic_dim)**2)
229 this%offset =
m_one/modq2
237 this%offset = this%offset*(-
expm1(-modq2/((
m_two*mu)**2)))
256 call pkernel_free(this%kernel)
257 call f_lib_finalize()
258 call dict_free(this%inputs)
267 class(
space_t),
intent(in) :: space
268 type(
cube_t),
intent(inout) :: cube
270 real(real64),
intent(in) :: qq_in(:)
272 real(real64) :: alpha, beta, gamma
273 real(real64) :: qq_abs(1:space%dim), qq(1:space%dim)
274 real(real64) :: modq2
287 call pkernel_free(this%kernel)
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)
296 call dict_set(this%inputs//
'kernel'//
'screening', coulb%mu)
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)
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)
305 call pkernel_set(this%kernel, verbose=
debug%info)
310 do idim = 1, space%periodic_dim
311 qq(idim) = qq_in(idim) - anint(qq_in(idim) +
m_half*1e-8_real64)
313 qq(space%periodic_dim + 1:space%dim) =
m_zero
315 modq2 = norm2(qq_abs)
316 if (modq2 > tol_vanishing_q)
then
319 this%offset = coulb%singularity
324 if (modq2 > tol_vanishing_q)
then
325 this%offset = this%offset*(-
expm1(-modq2/((
m_two*coulb%mu)**2)))
338 type(
mesh_t),
intent(in) :: mesh
340 real(real64),
intent(out) :: pot(:)
341 real(real64),
intent(in) :: rho(:)
347 type(coulomb_operator),
pointer :: kernel_pointer
348 real(real64) :: hartree_energy
349 real(real64) :: offset
356 real(real64),
allocatable :: pot_ion(:,:,:)
357 character(len=3) :: quiet
363 safe_allocate(pot_ion(1:cube%rs_n(1),1:cube%rs_n(2),1:cube%rs_n(3)))
366 if (
debug%info) quiet =
"NO"
380 kernel_pointer => this%kernel
381 call h_potential(this%datacode, kernel_pointer, cf%dRS, pot_ion, hartree_energy, offset, .false., quiet = quiet)
384 safe_deallocate_a(pot_ion)
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
402 character(len=3) :: quiet
406 type(coulomb_operator),
pointer :: kernel_pointer
407 real(real64) :: hartree_energy
408 real(real64) :: offset
416 real(real64),
allocatable :: pot_ion(:,:,:)
422 if (
present(sm))
then
428 safe_allocate(pot_ion(1:cube%rs_n(1), 1:cube%rs_n(2), 1:cube%rs_n(3)))
431 if (
debug%info) quiet =
"NO"
446 kernel_pointer => this%kernel
447 call h_potential(this%datacode, kernel_pointer, &
448 cf%dRS, pot_ion, hartree_energy, offset, .false., &
451 safe_deallocate_a(pot_ion)
453 if (
present(sm))
then
466 type(
cube_t),
intent(inout) :: cube
500 logical :: use_gradient
502 logical :: use_wb_corr
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 /)
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)
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)
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
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
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
subroutine, public kpoints_to_absolute(latt, kin, kout)
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 messages_not_implemented(feature, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_experimental(name, namespace)
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.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Describes mesh distribution to nodes.
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...