Octopus
poisson.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2011 M. Marques, A. Castro, A. Rubio,
2!! G. Bertsch, M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module poisson_oct_m
23 use accel_oct_m
24 use batch_oct_m
27 use comm_oct_m
28 use cube_oct_m
30 use debug_oct_m
32 use fft_oct_m
34 use global_oct_m
35 use index_oct_m
36 use, intrinsic :: iso_fortran_env
39 use math_oct_m
40 use mesh_oct_m
44 use mpi_oct_m
47#ifdef HAVE_OPENMP
48 use omp_lib
49#endif
51 use parser_oct_m
62 use space_oct_m
65 use types_oct_m
68 use xc_cam_oct_m
69
70 implicit none
71
72 private
73 public :: &
74 poisson_t, &
93
94 integer, public, parameter :: &
95 POISSON_DIRECT_SUM = -1, &
96 poisson_fft = 0, &
97 poisson_cg = 5, &
100 poisson_isf = 8, &
101 poisson_psolver = 10, &
102 poisson_no = -99, &
103 poisson_null = -999
104
105 type poisson_t
106 private
107 type(derivatives_t), pointer, public :: der
108 integer, public :: method = poisson_null
109 integer, public :: kernel
110 type(cube_t), public :: cube
111 type(mesh_cube_parallel_map_t), public :: mesh_cube_map
112 type(poisson_mg_solver_t) :: mg
113 type(poisson_fft_t), public :: fft_solver
114 real(real64), public :: poisson_soft_coulomb_param
115 logical :: all_nodes_default
116 type(poisson_corr_t) :: corrector
117 type(poisson_isf_t) :: isf_solver
118 type(poisson_psolver_t) :: psolver_solver
119 type(poisson_no_t) :: no_solver
120 integer :: nslaves
121 logical, public :: is_dressed = .false.
122 type(photon_mode_t), public :: photons
123#ifdef HAVE_MPI
124 type(MPI_Comm) :: intercomm
125 type(mpi_grp_t) :: local_grp
126 logical :: root
127#endif
128 end type poisson_t
129
130 integer, parameter :: &
131 CMD_FINISH = 1, &
133
134contains
135
136 !-----------------------------------------------------------------
137 subroutine poisson_init(this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx)
138 type(poisson_t), intent(inout) :: this
139 class(space_t), intent(in) :: space
140 type(namespace_t), intent(in) :: namespace
141 type(derivatives_t), target, intent(in) :: der
142 type(multicomm_t), intent(in) :: mc
143 type(stencil_t), intent(in) :: stencil
144 real(real64), optional, intent(in) :: qtot
145 character(len=*), optional, intent(in) :: label
146 integer, optional, intent(in) :: solver
147 logical, optional, intent(in) :: verbose
148 logical, optional, intent(in) :: force_serial
149 logical, optional, intent(in) :: force_cmplx
150
151 logical :: need_cube, isf_data_is_parallel
152 integer :: default_solver, default_kernel, box(space%dim), fft_type, fft_library
153 real(real64) :: fft_alpha
154 character(len=60) :: str
155
156 ! Make sure we do not try to initialize an already initialized solver
157 assert(this%method == poisson_null)
158
159 push_sub(poisson_init)
160
161 if (optional_default(verbose,.true.)) then
162 str = "Hartree"
163 if (present(label)) str = trim(label)
164 call messages_print_with_emphasis(msg=trim(str), namespace=namespace)
165 end if
166
167 this%nslaves = 0
168 this%der => der
169
170 !%Variable DressedOrbitals
171 !%Type logical
172 !%Default false
173 !%Section Hamiltonian::Poisson
174 !%Description
175 !% Allows for the calculation of coupled elecron-photon problems
176 !% by applying the dressed orbital approach. Details can be found in
177 !% https://arxiv.org/abs/1812.05562
178 !% At the moment, N electrons in d (<=3) spatial dimensions, coupled
179 !% to one photon mode can be described. The photon mode is included by
180 !% raising the orbital dimension to d+1 and changing the particle interaction
181 !% kernel and the local potential, where the former is included automatically,
182 !% but the latter needs to by added by hand as a user_defined_potential!
183 !% Coordinate 1-d: electron; coordinate d+1: photon.
184 !%End
185 call parse_variable(namespace, 'DressedOrbitals', .false., this%is_dressed)
186 call messages_print_var_value('DressedOrbitals', this%is_dressed, namespace=namespace)
187 if (this%is_dressed) then
188 assert(present(qtot))
189 call messages_experimental('Dressed Orbitals', namespace=namespace)
190 assert(qtot > m_zero)
191 call photon_mode_init(this%photons, namespace, der%dim-1)
192 call photon_mode_set_n_electrons(this%photons, qtot)
193 if(.not.allocated(this%photons%pol_dipole)) then
194 call photon_mode_compute_dipoles(this%photons, der%mesh)
195 end if
196 if (this%photons%nmodes > 1) then
197 call messages_not_implemented('DressedOrbitals for more than one photon mode', namespace=namespace)
198 end if
199 end if
201 this%all_nodes_default = .false.
202#ifdef HAVE_MPI
203 if (.not. optional_default(force_serial, .false.)) then
204 !%Variable ParallelizationPoissonAllNodes
205 !%Type logical
206 !%Default true
207 !%Section Execution::Parallelization
208 !%Description
209 !% When running in parallel, this variable selects whether the
210 !% Poisson solver should divide the work among all nodes or only
211 !% among the parallelization-in-domains groups.
212 !%End
214 call parse_variable(namespace, 'ParallelizationPoissonAllNodes', .true., this%all_nodes_default)
215 end if
216#endif
218 !%Variable PoissonSolver
219 !%Type integer
220 !%Section Hamiltonian::Poisson
221 !%Description
222 !% Defines which method to use to solve the Poisson equation. Some incompatibilities apply depending on
223 !% dimensionality, periodicity, etc.
224 !% For a comparison of the accuracy and performance of the methods in Octopus, see P Garcia-Risue&ntilde;o,
225 !% J Alberdi-Rodriguez <i>et al.</i>, <i>J. Comp. Chem.</i> <b>35</b>, 427-444 (2014)
226 !% or <a href=http://arxiv.org/abs/1211.2092>arXiV</a>.
227 !% Defaults:
228 !% <br> 1D and 2D: <tt>fft</tt>.
229 !% <br> 3D: <tt>cg_corrected</tt> if curvilinear, <tt>isf</tt> if not periodic, <tt>fft</tt> if periodic.
230 !% <br> Dressed orbitals: <tt>direct_sum</tt>.
231 !%Option NoPoisson -99
232 !% Do not use a Poisson solver at all.
233 !%Option direct_sum -1
234 !% Direct evaluation of the Hartree potential (only for finite systems).
235 !%Option fft 0
236 !% The Poisson equation is solved using FFTs. A cutoff technique
237 !% for the Poisson kernel is selected so the proper boundary
238 !% conditions are imposed according to the periodicity of the
239 !% system. This can be overridden by the <tt>PoissonFFTKernel</tt>
240 !% variable. To choose the FFT library use <tt>FFTLibrary</tt>
241 !%Option cg 5
242 !% Conjugate gradients (only for finite systems).
243 !%Option cg_corrected 6
244 !% Conjugate gradients, corrected for boundary conditions (only for finite systems).
245 !%Option multigrid 7
246 !% Multigrid method (only for finite systems).
247 !%Option isf 8
248 !% Interpolating Scaling Functions Poisson solver (only for finite systems).
249 !%Option psolver 10
250 !% Solver based on Interpolating Scaling Functions as implemented in the PSolver library.
251 !% Parallelization in k-points requires <tt>PoissonSolverPSolverParallelData</tt> = no.
252 !% Requires the PSolver external library.
253 !%End
254
255 default_solver = poisson_fft
256
257 if (space%dim == 3 .and. .not. space%is_periodic()) default_solver = poisson_isf
258
259#ifdef HAVE_CUDA
260 if(accel_is_enabled()) default_solver = poisson_fft
261#endif
262
263 if (space%dim > 3) default_solver = poisson_no ! Kernel for higher dimensions is not implemented.
264
265 if (der%mesh%use_curvilinear) then
266 select case (space%dim)
267 case (1)
268 default_solver = poisson_direct_sum
269 case (2)
270 default_solver = poisson_direct_sum
271 case (3)
272 default_solver = poisson_multigrid
273 end select
274 end if
275
276 if (this%is_dressed) default_solver = poisson_direct_sum
277
278 if (.not. present(solver)) then
279 call parse_variable(namespace, 'PoissonSolver', default_solver, this%method)
280 else
281 this%method = solver
282 end if
283 if (.not. varinfo_valid_option('PoissonSolver', this%method)) call messages_input_error(namespace, 'PoissonSolver')
284 if (optional_default(verbose, .true.)) then
285 select case (this%method)
286 case (poisson_direct_sum)
287 str = "direct sum"
288 case (poisson_fft)
289 str = "fast Fourier transform"
290 case (poisson_cg)
291 str = "conjugate gradients"
293 str = "conjugate gradients, corrected"
294 case (poisson_multigrid)
295 str = "multigrid"
296 case (poisson_isf)
297 str = "interpolating scaling functions"
298 case (poisson_psolver)
299 str = "interpolating scaling functions (from BigDFT)"
300 case (poisson_no)
301 str = "no Poisson solver - Hartree set to 0"
302 end select
303 write(message(1),'(a,a,a)') "The chosen Poisson solver is '", trim(str), "'"
304 call messages_info(1, namespace=namespace)
305 end if
306
307 if (space%dim > 3 .and. this%method /= poisson_no) then
308 call messages_input_error(namespace, 'PoissonSolver', 'Currently no Poisson solver is available for Dimensions > 3')
309 end if
310
311 if (this%method /= poisson_fft) then
312 this%kernel = poisson_fft_kernel_none
313 else
314
315 ! Documentation in cube.F90
316 call parse_variable(namespace, 'FFTLibrary', fftlib_fftw, fft_library)
317
318 !%Variable PoissonFFTKernel
319 !%Type integer
320 !%Section Hamiltonian::Poisson
321 !%Description
322 !% Defines which kernel is used to impose the correct boundary
323 !% conditions when using FFTs to solve the Poisson equation. The
324 !% default is selected depending on the dimensionality and
325 !% periodicity of the system:
326 !% <br>In 1D, <tt>spherical</tt> if finite, <tt>fft_nocut</tt> if periodic.
327 !% <br>In 2D, <tt>spherical</tt> if finite, <tt>cylindrical</tt> if 1D-periodic, <tt>fft_nocut</tt> if 2D-periodic.
328 !% <br>In 3D, <tt>spherical</tt> if finite, <tt>cylindrical</tt> if 1D-periodic, <tt>planar</tt> if 2D-periodic,
329 !% <tt>fft_nocut</tt> if 3D-periodic.
330 !% See C. A. Rozzi et al., <i>Phys. Rev. B</i> <b>73</b>, 205119 (2006) for 3D implementation and
331 !% A. Castro et al., <i>Phys. Rev. B</i> <b>80</b>, 033102 (2009) for 2D implementation.
332 !%Option spherical 0
333 !% FFTs using spherical cutoff (in 2D or 3D).
334 !%Option cylindrical 1
335 !% FFTs using cylindrical cutoff (in 2D or 3D).
336 !%Option planar 2
337 !% FFTs using planar cutoff (in 3D).
338 !%Option fft_nocut 3
339 !% FFTs without using a cutoff (for fully periodic systems).
340 !%Option multipole_correction 4
341 !% The boundary conditions are imposed by using a multipole expansion. Only appropriate for finite systems.
342 !% Further specification occurs with variables <tt>PoissonSolverBoundaries</tt> and <tt>PoissonSolverMaxMultipole</tt>.
343 !%End
344
345 select case (space%dim)
346 case (1)
347 if (.not. space%is_periodic()) then
348 default_kernel = poisson_fft_kernel_sph
349 else
350 default_kernel = poisson_fft_kernel_nocut
351 end if
352 case (2)
353 if (space%periodic_dim == 2) then
354 default_kernel = poisson_fft_kernel_nocut
355 else if (space%is_periodic()) then
356 default_kernel = space%periodic_dim
357 else
358 default_kernel = poisson_fft_kernel_sph
359 end if
360 case (3)
361 default_kernel = space%periodic_dim
362 end select
363
364 call parse_variable(namespace, 'PoissonFFTKernel', default_kernel, this%kernel)
365 if (.not. varinfo_valid_option('PoissonFFTKernel', this%kernel)) call messages_input_error(namespace, 'PoissonFFTKernel')
366
367 if (optional_default(verbose,.true.)) then
368 call messages_print_var_option("PoissonFFTKernel", this%kernel, namespace=namespace)
369 end if
370
371 ! the multipole correction kernel does not work on GPUs
372 if(this%kernel == poisson_fft_kernel_corrected .and. fft_default_lib == fftlib_accel) then
374 message(1) = 'PoissonFFTKernel=multipole_correction is not supported on GPUs'
375 message(2) = 'Using FFTW to compute the FFTs on the CPU'
376 call messages_info(2, namespace=namespace)
377 end if
378
379 end if
380
381 !We assume the developer knows what he is doing by providing the solver option
382 if (.not. present(solver)) then
383 if (space%is_periodic() .and. this%method == poisson_direct_sum) then
384 message(1) = 'A periodic system may not use the direct_sum Poisson solver.'
385 call messages_fatal(1, namespace=namespace)
386 end if
387
388 if (space%is_periodic() .and. this%method == poisson_cg_corrected) then
389 message(1) = 'A periodic system may not use the cg_corrected Poisson solver.'
390 call messages_fatal(1, namespace=namespace)
391 end if
392
393 if (space%is_periodic() .and. this%method == poisson_cg) then
394 message(1) = 'A periodic system may not use the cg Poisson solver.'
395 call messages_fatal(1, namespace=namespace)
396 end if
397
398 if (space%is_periodic() .and. this%method == poisson_multigrid) then
399 message(1) = 'A periodic system may not use the multigrid Poisson solver.'
400 call messages_fatal(1, namespace=namespace)
401 end if
402
403 select case (space%dim)
404 case (1)
405
406 select case (space%periodic_dim)
407 case (0)
408 if ((this%method /= poisson_fft) .and. (this%method /= poisson_direct_sum)) then
409 message(1) = 'A finite 1D system may only use fft or direct_sum Poisson solvers.'
410 call messages_fatal(1, namespace=namespace)
411 end if
412 case (1)
413 if (this%method /= poisson_fft) then
414 message(1) = 'A periodic 1D system may only use the fft Poisson solver.'
415 call messages_fatal(1, namespace=namespace)
416 end if
417 end select
418
419 if (der%mesh%use_curvilinear .and. this%method /= poisson_direct_sum) then
420 message(1) = 'If curvilinear coordinates are used in 1D, then the only working'
421 message(2) = 'Poisson solver is direct_sum.'
422 call messages_fatal(2, namespace=namespace)
423 end if
424
425 case (2)
426
427 if ((this%method /= poisson_fft) .and. (this%method /= poisson_direct_sum)) then
428 message(1) = 'A 2D system may only use fft or direct_sum solvers.'
429 call messages_fatal(1, namespace=namespace)
430 end if
431
432 if (der%mesh%use_curvilinear .and. (this%method /= poisson_direct_sum)) then
433 message(1) = 'If curvilinear coordinates are used in 2D, then the only working'
434 message(2) = 'Poisson solver is direct_sum.'
435 call messages_fatal(2, namespace=namespace)
436 end if
437
438 case (3)
439
440 if (space%is_periodic() .and. this%method == poisson_isf) then
441 call messages_write('The ISF solver can only be used for finite systems.')
442 call messages_fatal()
443 end if
444
445 if (space%is_periodic() .and. this%method == poisson_fft .and. &
446 this%kernel /= space%periodic_dim .and. this%kernel >= 0 .and. this%kernel <= 3) then
447 write(message(1), '(a,i1,a)')'The system is periodic in ', space%periodic_dim ,' dimension(s),'
448 write(message(2), '(a,i1,a)')'but Poisson solver is set for ', this%kernel, ' dimensions.'
449 call messages_warning(2, namespace=namespace)
450 end if
451
452 if (space%is_periodic() .and. this%method == poisson_fft .and. this%kernel == poisson_fft_kernel_corrected) then
453 write(message(1), '(a,i1,a)')'PoissonFFTKernel = multipole_correction cannot be used for periodic systems.'
454 call messages_fatal(1, namespace=namespace)
455 end if
456
457 if (der%mesh%use_curvilinear .and. .not. any(this%method == [poisson_cg_corrected, poisson_multigrid])) then
458 message(1) = 'If curvilinear coordinates are used, then the only working'
459 message(2) = 'Poisson solver are cg_corrected and multigrid.'
460 call messages_fatal(2, namespace=namespace)
461 end if
462 if (der%mesh%use_curvilinear .and. this%method == poisson_multigrid .and. accel_is_enabled()) then
463 call messages_not_implemented('Multigrid Poisson solver with curvilinear coordinates on GPUs')
464 end if
465
466 select type (box => der%mesh%box)
467 type is (box_minimum_t)
468 if (this%method == poisson_cg_corrected) then
469 message(1) = 'When using the "minimum" box shape and the "cg_corrected"'
470 message(2) = 'Poisson solver, we have observed "sometimes" some non-'
471 message(3) = 'negligible error. You may want to check that the "fft" or "cg"'
472 message(4) = 'solver are providing, in your case, the same results.'
473 call messages_warning(4, namespace=namespace)
474 end if
475 end select
476
477 end select
478 end if
479
480 if (this%method == poisson_psolver) then
481#if !(defined HAVE_PSOLVER)
482 message(1) = "The PSolver Poisson solver cannot be used since the code was not compiled with the PSolver library."
483 call messages_fatal(1, namespace=namespace)
484#endif
485 end if
486
487 if (optional_default(verbose,.true.)) then
488 call messages_print_with_emphasis(namespace=namespace)
489 end if
490
491 ! Now that we know the method, we check if we need a cube and its dimentions
492 need_cube = .false.
493 fft_type = fft_real
494 if (optional_default(force_cmplx, .false.)) fft_type = fft_complex
495
496 if (this%method == poisson_isf .or. this%method == poisson_psolver) then
497 fft_type = fft_none
498 box(:) = der%mesh%idx%ll(:)
499 need_cube = .true.
500 end if
501
502 if (this%method == poisson_psolver .and. multicomm_have_slaves(mc)) then
503 call messages_not_implemented('Task parallelization with PSolver Poisson solver', namespace=namespace)
504 end if
505
507 ! Documentation in poisson_psolver.F90
508 call parse_variable(namespace, 'PoissonSolverPSolverParallelData', .true., isf_data_is_parallel)
509 if (this%method == poisson_psolver .and. isf_data_is_parallel) then
510 call messages_not_implemented("k-point parallelization with PSolver library and", namespace=namespace)
511 call messages_not_implemented("PoissonSolverPSolverParallelData = yes", namespace=namespace)
512 end if
513 if (this%method == poisson_fft .and. fft_library == fftlib_pfft) then
514 call messages_not_implemented("k-point parallelization with PFFT library for", namespace=namespace)
515 call messages_not_implemented("PFFT library for Poisson solver", namespace=namespace)
516 end if
517 end if
518
519 if (this%method == poisson_fft) then
520
521 need_cube = .true.
522
523 !%Variable DoubleFFTParameter
524 !%Type float
525 !%Default 2.0
526 !%Section Mesh::FFTs
527 !%Description
528 !% For solving the Poisson equation in Fourier space, and for applying the local potential
529 !% in Fourier space, an auxiliary cubic mesh is built. This mesh will be larger than
530 !% the circumscribed cube of the usual mesh by a factor <tt>DoubleFFTParameter</tt>. See
531 !% the section that refers to Poisson equation, and to the local potential for details
532 !% [the default value of two is typically good].
533 !%End
534 call parse_variable(namespace, 'DoubleFFTParameter', m_two, fft_alpha)
535 if (fft_alpha < m_one .or. fft_alpha > m_three) then
536 write(message(1), '(a,f12.5,a)') "Input: '", fft_alpha, &
537 "' is not a valid DoubleFFTParameter"
538 message(2) = '1.0 <= DoubleFFTParameter <= 3.0'
539 call messages_fatal(2, namespace=namespace)
540 end if
541
542 if (space%dim /= 3 .and. fft_library == fftlib_pfft) then
543 call messages_not_implemented('PFFT support for dimensionality other than 3', namespace=namespace)
544 end if
545
546 select case (space%dim)
547
548 case (1)
549 select case (this%kernel)
551 call mesh_double_box(space, der%mesh, fft_alpha, box)
553 box = der%mesh%idx%ll
554 end select
555
556 case (2)
557 select case (this%kernel)
559 call mesh_double_box(space, der%mesh, fft_alpha, box)
560 box(1:2) = maxval(box)
562 call mesh_double_box(space, der%mesh, fft_alpha, box)
564 box(:) = der%mesh%idx%ll(:)
565 end select
566
567 case (3)
568 select case (this%kernel)
570 call mesh_double_box(space, der%mesh, fft_alpha, box)
571 box(:) = maxval(box)
573 call mesh_double_box(space, der%mesh, fft_alpha, box)
574 box(2) = maxval(box(2:3)) ! max of finite directions
575 box(3) = maxval(box(2:3)) ! max of finite directions
577 box(:) = der%mesh%idx%ll(:)
579 call mesh_double_box(space, der%mesh, fft_alpha, box)
580 end select
581
582 end select
583
584 end if
585
586 ! Create the cube
587 if (need_cube) then
588 call cube_init(this%cube, box, namespace, space, der%mesh%spacing, &
589 der%mesh%coord_system, fft_type = fft_type, &
590 need_partition=.not.der%mesh%parallel_in_domains)
591 call cube_init_cube_map(this%cube, der%mesh)
592 if (this%cube%parallel_in_domains .and. this%method == poisson_fft) then
593 call mesh_cube_parallel_map_init(this%mesh_cube_map, der%mesh, this%cube)
594 end if
595 end if
596
597 if (this%is_dressed .and. .not. this%method == poisson_direct_sum) then
598 write(message(1), '(a)')'Dressed Orbital calculation currently only working with direct sum Poisson solver.'
599 call messages_fatal(1, namespace=namespace)
600 end if
601
602 call poisson_kernel_init(this, namespace, space, mc, stencil)
603
604 pop_sub(poisson_init)
605 end subroutine poisson_init
606
607 !-----------------------------------------------------------------
608 subroutine poisson_end(this)
609 type(poisson_t), intent(inout) :: this
610
611 logical :: has_cube
612
613 push_sub(poisson_end)
614
615 has_cube = .false.
616
617 select case (this%method)
618 case (poisson_fft)
619 call poisson_fft_end(this%fft_solver)
620 if (this%kernel == poisson_fft_kernel_corrected) call poisson_corrections_end(this%corrector)
621 has_cube = .true.
622
624 call poisson_cg_end()
625 call poisson_corrections_end(this%corrector)
626
627 case (poisson_multigrid)
628 call poisson_multigrid_end(this%mg)
629
630 case (poisson_isf)
631 call poisson_isf_end(this%isf_solver)
632 has_cube = .true.
633
634 case (poisson_psolver)
635 call poisson_psolver_end(this%psolver_solver)
636 has_cube = .true.
637
638 case (poisson_no)
639 call poisson_no_end(this%no_solver)
640
641 end select
642 this%method = poisson_null
643
644 if (has_cube) then
645 if (this%cube%parallel_in_domains) then
646 call mesh_cube_parallel_map_end(this%mesh_cube_map)
647 end if
648 call cube_end(this%cube)
649 end if
650
651 if (this%is_dressed) then
652 call photon_mode_end(this%photons)
653 end if
654 this%is_dressed = .false.
655
656 pop_sub(poisson_end)
657 end subroutine poisson_end
658
659 !-----------------------------------------------------------------
660
661 subroutine zpoisson_solve_real_and_imag_separately(this, namespace, pot, rho, all_nodes, kernel)
662 type(poisson_t), intent(in) :: this
663 type(namespace_t), intent(in) :: namespace
664 complex(real64), contiguous, intent(inout) :: pot(:)
665 complex(real64), contiguous, intent(in) :: rho(:)
666 logical, optional, intent(in) :: all_nodes
667 type(fourier_space_op_t), optional, intent(in) :: kernel
668
669 real(real64), allocatable :: aux1(:), aux2(:)
670 type(derivatives_t), pointer :: der
671 logical :: all_nodes_value
672
673
674 der => this%der
675
677
678 call profiling_in('POISSON_RE_IM_SOLVE')
679
680 if (present(kernel)) then
681 assert(.not. any(abs(kernel%qq(:))>1e-8_real64))
682 end if
683
684 all_nodes_value = optional_default(all_nodes, this%all_nodes_default)
685
686 safe_allocate(aux1(1:der%mesh%np))
687 safe_allocate(aux2(1:der%mesh%np))
688 ! first the real part
689 aux1(1:der%mesh%np) = real(rho(1:der%mesh%np), real64)
690 aux2(1:der%mesh%np) = real(pot(1:der%mesh%np), real64)
691 call dpoisson_solve(this, namespace, aux2, aux1, all_nodes=all_nodes_value, kernel=kernel)
692 pot(1:der%mesh%np) = aux2(1:der%mesh%np)
693
694 ! now the imaginary part
695 aux1(1:der%mesh%np) = aimag(rho(1:der%mesh%np))
696 aux2(1:der%mesh%np) = aimag(pot(1:der%mesh%np))
697 call dpoisson_solve(this, namespace, aux2, aux1, all_nodes=all_nodes_value, kernel=kernel)
698 pot(1:der%mesh%np) = pot(1:der%mesh%np) + m_zi*aux2(1:der%mesh%np)
699
700 safe_deallocate_a(aux1)
701 safe_deallocate_a(aux2)
702
703 call profiling_out('POISSON_RE_IM_SOLVE')
704
707
708 !-----------------------------------------------------------------
709
710 subroutine zpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
711 type(poisson_t), intent(in) :: this
712 type(namespace_t), intent(in) :: namespace
713 complex(real64), contiguous, intent(inout) :: pot(:)
714 complex(real64), contiguous, intent(in) :: rho(:)
715 logical, optional, intent(in) :: all_nodes
716 type(fourier_space_op_t), optional, intent(in) :: kernel
717 logical, optional, intent(in) :: reset
718
719 logical :: all_nodes_value
720
721 push_sub(zpoisson_solve)
722
723 all_nodes_value = optional_default(all_nodes, this%all_nodes_default)
724
725 assert(ubound(pot, dim = 1) == this%der%mesh%np_part .or. ubound(pot, dim = 1) == this%der%mesh%np)
726 assert(ubound(rho, dim = 1) == this%der%mesh%np_part .or. ubound(rho, dim = 1) == this%der%mesh%np)
727
728 assert(this%method /= poisson_null)
729
730 if (poisson_solver_is_iterative(this) .and. optional_default(reset, .true.)) then
731 pot(1:this%der%mesh%np) = m_zero
732 end if
733
734 if (this%method == poisson_fft .and. this%kernel /= poisson_fft_kernel_corrected &
735 .and. .not. this%is_dressed) then
736 !The default (real) Poisson solver is used for OEP and Sternheimer calls were we do not need
737 !a complex-to-xomplex FFT as these parts use the normal Coulomb potential
738 if (this%cube%fft%type == fft_complex) then
739 !We add the profiling here, as the other path uses dpoisson_solve
740 call profiling_in('ZPOISSON_SOLVE')
741 call zpoisson_fft_solve(this%fft_solver, this%der%mesh, this%cube, pot, rho, this%mesh_cube_map, kernel=kernel)
742 call profiling_out('ZPOISSON_SOLVE')
743 else
744 call zpoisson_solve_real_and_imag_separately(this, namespace, pot, rho, all_nodes_value, kernel=kernel)
745 end if
746 else
747 call zpoisson_solve_real_and_imag_separately(this, namespace, pot, rho, all_nodes_value, kernel = kernel)
748 end if
749
750 pop_sub(zpoisson_solve)
751 end subroutine zpoisson_solve
752
753
754 !-----------------------------------------------------------------
755
756 subroutine poisson_solve_batch(this, namespace, potb, rhob, all_nodes, kernel)
757 type(poisson_t), intent(inout) :: this
758 type(namespace_t), intent(in) :: namespace
759 type(batch_t), intent(inout) :: potb
760 type(batch_t), intent(inout) :: rhob
761 logical, optional, intent(in) :: all_nodes
762 type(fourier_space_op_t), optional, intent(in) :: kernel
763
764 integer :: ii
765
766 push_sub(poisson_solve_batch)
767
768 assert(potb%nst_linear == rhob%nst_linear)
769 assert(potb%type() == rhob%type())
770
771 if (potb%type() == type_float) then
772 do ii = 1, potb%nst_linear
773 call dpoisson_solve(this, namespace, potb%dff_linear(:, ii), rhob%dff_linear(:, ii), all_nodes, kernel=kernel)
774 end do
775 else
776 do ii = 1, potb%nst_linear
777 call zpoisson_solve(this, namespace, potb%zff_linear(:, ii), rhob%zff_linear(:, ii), all_nodes, kernel=kernel)
778 end do
779 end if
780
781 pop_sub(poisson_solve_batch)
782 end subroutine poisson_solve_batch
783
784 !-----------------------------------------------------------------
785
791 subroutine dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
792 type(poisson_t), intent(in) :: this
793 type(namespace_t), intent(in) :: namespace
794 real(real64), contiguous, intent(inout) :: pot(:)
795 real(real64), contiguous, intent(in) :: rho(:)
799 logical, optional, intent(in) :: all_nodes
800 type(fourier_space_op_t), optional, intent(in) :: kernel
801 logical, optional, intent(in) :: reset
802
803 type(derivatives_t), pointer :: der
804 real(real64), allocatable :: rho_corrected(:), vh_correction(:)
805 logical :: all_nodes_value
806
807 call profiling_in('POISSON_SOLVE')
808 push_sub(dpoisson_solve)
809
810 der => this%der
811
812 assert(ubound(pot, dim = 1) == der%mesh%np_part .or. ubound(pot, dim = 1) == der%mesh%np)
813 assert(ubound(rho, dim = 1) == der%mesh%np_part .or. ubound(rho, dim = 1) == der%mesh%np)
814
815 ! Check optional argument and set to default if necessary.
816 all_nodes_value = optional_default(all_nodes, this%all_nodes_default)
817
818 if (poisson_solver_is_iterative(this) .and. optional_default(reset, .true.)) then
819 pot(1:der%mesh%np) = m_zero
820 end if
821
822 assert(this%method /= poisson_null)
823
824 if (present(kernel)) then
825 assert(this%method == poisson_fft)
826 end if
827
828 select case (this%method)
829 case (poisson_direct_sum)
830 if ((this%is_dressed .and. this%der%dim - 1 > 3) .or. this%der%dim > 3) then
831 message(1) = "Direct sum Poisson solver only available for 1, 2, or 3 dimensions."
832 call messages_fatal(1, namespace=namespace)
833 end if
834 call poisson_solve_direct(this, namespace, pot, rho)
835
837 call poisson_cg1(namespace, der, this%corrector, pot, rho)
838
840 safe_allocate(rho_corrected(1:der%mesh%np))
841 safe_allocate(vh_correction(1:der%mesh%np_part))
842
843 call correct_rho(this%corrector, der, rho, rho_corrected, vh_correction)
844
845 pot(1:der%mesh%np) = pot(1:der%mesh%np) - vh_correction(1:der%mesh%np)
846 call poisson_cg2(namespace, der, pot, rho_corrected)
847 pot(1:der%mesh%np) = pot(1:der%mesh%np) + vh_correction(1:der%mesh%np)
848
849 safe_deallocate_a(rho_corrected)
850 safe_deallocate_a(vh_correction)
851
852 case (poisson_multigrid)
853 call poisson_multigrid_solver(this%mg, namespace, der, pot, rho)
854
855 case (poisson_fft)
856 if (this%kernel /= poisson_fft_kernel_corrected) then
857 call dpoisson_fft_solve(this%fft_solver, der%mesh, this%cube, pot, rho, this%mesh_cube_map, kernel=kernel)
858 else
859 safe_allocate(rho_corrected(1:der%mesh%np))
860 safe_allocate(vh_correction(1:der%mesh%np_part))
861
862 call correct_rho(this%corrector, der, rho, rho_corrected, vh_correction)
863 call dpoisson_fft_solve(this%fft_solver, der%mesh, this%cube, pot, rho_corrected, this%mesh_cube_map, &
864 average_to_zero = .true., kernel=kernel)
865
866 pot(1:der%mesh%np) = pot(1:der%mesh%np) + vh_correction(1:der%mesh%np)
867 safe_deallocate_a(rho_corrected)
868 safe_deallocate_a(vh_correction)
869 end if
870
872 call poisson_isf_solve(this%isf_solver, der%mesh, this%cube, pot, rho, all_nodes_value)
873
874
875 case (poisson_psolver)
876 if (this%psolver_solver%datacode == "G") then
877 ! Global version
878 call poisson_psolver_global_solve(this%psolver_solver, der%mesh, this%cube, pot, rho)
879 else ! "D" Distributed version
880 call poisson_psolver_parallel_solve(this%psolver_solver, der%mesh, this%cube, pot, rho, this%mesh_cube_map)
881 end if
882
883 case (poisson_no)
884 call poisson_no_solve(this%no_solver, der%mesh, pot, rho)
885 end select
886
887
888 ! Add extra terms for dressed interaction
889 if (this%is_dressed .and. this%method /= poisson_no) then
890 call photon_mode_add_poisson_terms(this%photons, der%mesh, rho, pot)
891 end if
892
893 pop_sub(dpoisson_solve)
894 call profiling_out('POISSON_SOLVE')
895 end subroutine dpoisson_solve
896
897 !-----------------------------------------------------------------
898 subroutine poisson_init_sm(this, namespace, space, main, der, sm, method, force_cmplx)
899 type(poisson_t), intent(inout) :: this
900 type(namespace_t), intent(in) :: namespace
901 class(space_t), intent(in) :: space
902 type(poisson_t), intent(in) :: main
903 type(derivatives_t), target, intent(in) :: der
904 type(submesh_t), intent(inout) :: sm
905 integer, optional, intent(in) :: method
906 logical, optional, intent(in) :: force_cmplx
907
908 integer :: default_solver, idir, iter, maxl
909 integer :: box(space%dim)
910 real(real64) :: qq(der%dim), threshold
911
912 if (this%method /= poisson_null) return ! already initialized
913
914 push_sub(poisson_init_sm)
915
916 this%is_dressed = .false.
917 !TODO: To be implemented as an option
918 this%all_nodes_default = .false.
919
920 this%nslaves = 0
921 this%der => der
922
923#ifdef HAVE_MPI
924 this%all_nodes_default = main%all_nodes_default
925#endif
926
927 default_solver = poisson_direct_sum
928 this%method = default_solver
929 if (present(method)) this%method = method
930
931 if (der%mesh%use_curvilinear) then
932 call messages_not_implemented("Submesh Poisson solver with curvilinear mesh", namespace=namespace)
933 end if
934
935 this%kernel = poisson_fft_kernel_none
936
937 select case (this%method)
938 case (poisson_direct_sum)
939 !Nothing to be done
940
941 case (poisson_isf)
942 !TODO: Add support for domain parrallelization
943 assert(.not. der%mesh%parallel_in_domains)
944 call submesh_get_cube_dim(sm, space, box)
945 call submesh_init_cube_map(sm, space)
946 call cube_init(this%cube, box, namespace, space, sm%mesh%spacing, sm%mesh%coord_system, &
947 fft_type = fft_none, need_partition=.not.der%mesh%parallel_in_domains)
948 call cube_init_cube_map(this%cube, sm%mesh)
949 call poisson_isf_init(this%isf_solver, namespace, der%mesh, this%cube, mpi_world%comm, init_world = this%all_nodes_default)
950
951 case (poisson_psolver)
952 !TODO: Add support for domain parrallelization
953 assert(.not. der%mesh%parallel_in_domains)
954 if (this%all_nodes_default) then
955 this%cube%mpi_grp = mpi_world
956 else
957 this%cube%mpi_grp = this%der%mesh%mpi_grp
958 end if
959 call submesh_get_cube_dim(sm, space, box)
960 call submesh_init_cube_map(sm, space)
961 call cube_init(this%cube, box, namespace, space, sm%mesh%spacing, sm%mesh%coord_system, &
962 fft_type = fft_none, need_partition=.not.der%mesh%parallel_in_domains)
963 call cube_init_cube_map(this%cube, sm%mesh)
964 qq = m_zero
965 call poisson_psolver_init(this%psolver_solver, namespace, space, this%cube, m_zero, qq, force_isolated=.true.)
966 call poisson_psolver_get_dims(this%psolver_solver, this%cube)
967 case (poisson_fft)
968 !Here we impose zero boundary conditions
969 this%kernel = poisson_fft_kernel_sph
970 !We need to parse this, in case this routine is called before poisson_init
971 call parse_variable(namespace, 'FFTLibrary', fftlib_fftw, fft_default_lib)
972
973 call submesh_get_cube_dim(sm, space, box)
974 call submesh_init_cube_map(sm, space)
975 !We double the size of the cell
976 !Maybe the factor of two should be controlled as a variable
977 do idir = 1, space%dim
978 box(idir) = (2 * (box(idir) - 1)) + 1
979 end do
980 if (optional_default(force_cmplx, .false.)) then
981 call cube_init(this%cube, box, namespace, space, sm%mesh%spacing, sm%mesh%coord_system, &
982 fft_type = fft_complex, need_partition=.not.der%mesh%parallel_in_domains)
983 else
984 call cube_init(this%cube, box, namespace, space, sm%mesh%spacing, sm%mesh%coord_system, &
985 fft_type = fft_real, need_partition=.not.der%mesh%parallel_in_domains)
986 end if
987 call poisson_fft_init(this%fft_solver, namespace, space, this%cube, this%kernel)
988 case (poisson_cg)
989 call parse_variable(namespace, 'PoissonSolverMaxMultipole', 4, maxl)
990 write(message(1),'(a,i2)')'Info: Boundary conditions fixed up to L =', maxl
991 call messages_info(1, namespace=namespace)
992 call parse_variable(namespace, 'PoissonSolverMaxIter', 500, iter)
993 call parse_variable(namespace, 'PoissonSolverThreshold', 1.0e-6_real64, threshold)
994 call poisson_corrections_init(this%corrector, namespace, space, maxl, this%der%mesh)
995 call poisson_cg_init(threshold, iter)
996 end select
997
998 pop_sub(poisson_init_sm)
999 end subroutine poisson_init_sm
1000
1001 ! -----------------------------------------------------------------
1002
1003 logical pure function poisson_solver_is_iterative(this) result(iterative)
1004 type(poisson_t), intent(in) :: this
1005
1006 iterative = this%method == poisson_cg .or. this%method == poisson_cg_corrected
1007 end function poisson_solver_is_iterative
1008
1009 !-----------------------------------------------------------------
1010 subroutine poisson_async_init(this, mc)
1011 type(poisson_t), intent(inout) :: this
1012 type(multicomm_t), intent(in) :: mc
1013
1014 push_sub(poisson_async_init)
1015
1016#ifdef HAVE_MPI
1017 if (multicomm_have_slaves(mc)) then
1018
1019 call mpi_grp_init(this%local_grp, mc%group_comm(p_strategy_states))
1020
1021 this%root = (this%local_grp%is_root())
1022
1023 this%intercomm = mc%slave_intercomm
1024 call mpi_comm_remote_size(this%intercomm, this%nslaves)
1025
1026 end if
1027#endif
1028
1029 pop_sub(poisson_async_init)
1030
1031 end subroutine poisson_async_init
1032
1033 !-----------------------------------------------------------------
1034
1035 subroutine poisson_async_end(this, mc)
1036 type(poisson_t), intent(inout) :: this
1037 type(multicomm_t), intent(in) :: mc
1038
1039#ifdef HAVE_MPI
1040 integer :: islave
1041#endif
1042
1043 push_sub(poisson_async_end)
1044
1045#ifdef HAVE_MPI
1046 if (multicomm_have_slaves(mc)) then
1047
1048 ! send the finish signal
1049 do islave = this%local_grp%rank, this%nslaves - 1, this%local_grp%size
1050 call mpi_send(m_one, 1, mpi_double_precision, islave, cmd_finish, this%intercomm)
1051 end do
1052
1053 end if
1054#endif
1055
1056 pop_sub(poisson_async_end)
1057
1058 end subroutine poisson_async_end
1059
1060 !-----------------------------------------------------------------
1061
1062 subroutine poisson_slave_work(this, namespace)
1063 type(poisson_t), intent(inout) :: this
1064 type(namespace_t), intent(in) :: namespace
1065
1066#ifdef HAVE_MPI
1067 real(real64), allocatable :: rho(:), pot(:)
1068 logical :: done
1069 type(mpi_status) :: status
1070 integer :: bcast_root
1071
1072 push_sub(poisson_slave_work)
1073 call profiling_in("SLAVE_WORK")
1074
1075 safe_allocate(rho(1:this%der%mesh%np))
1076 safe_allocate(pot(1:this%der%mesh%np))
1077 done = .false.
1078
1079 do while(.not. done)
1080
1081 call profiling_in("SLAVE_WAIT")
1082 call mpi_recv(rho(1), this%der%mesh%np, mpi_double_precision, mpi_any_source, mpi_any_tag, this%intercomm, status)
1083 call profiling_out("SLAVE_WAIT")
1084
1085 ! The tag of the message tells us what we have to do.
1086 select case (status%MPI_TAG)
1087
1088 case (cmd_finish)
1089 done = .true.
1091 case (cmd_poisson_solve)
1092 call dpoisson_solve(this, namespace, pot, rho)
1093
1094 call profiling_in("SLAVE_BROADCAST")
1095 bcast_root = mpi_proc_null
1096 if (this%root) bcast_root = mpi_root
1097 call mpi_bcast(pot(1), this%der%mesh%np, mpi_double_precision, bcast_root, this%intercomm)
1098 call profiling_out("SLAVE_BROADCAST")
1099
1100 end select
1101
1102 end do
1103
1104 safe_deallocate_a(pot)
1105 safe_deallocate_a(rho)
1106
1107 call profiling_out("SLAVE_WORK")
1108 pop_sub(poisson_slave_work)
1109#endif
1110 end subroutine poisson_slave_work
1111
1112 !----------------------------------------------------------------
1113
1114 logical pure function poisson_is_async(this) result(async)
1115 type(poisson_t), intent(in) :: this
1116
1117 async = (this%nslaves > 0)
1119 end function poisson_is_async
1120
1121 !----------------------------------------------------------------
1122
1123 subroutine poisson_build_kernel(this, namespace, space, coulb, qq, cam, singul)
1124 type(poisson_t), intent(in) :: this
1125 type(namespace_t), intent(in) :: namespace
1126 class(space_t), intent(in) :: space
1127 type(fourier_space_op_t), intent(inout) :: coulb
1128 real(real64), intent(in) :: qq(:)
1129 type(xc_cam_t), intent(in) :: cam
1130 real(real64), optional, intent(in) :: singul
1131
1132 real(real64), parameter :: vanishing_q = 1.0e-5_real64
1133 logical :: reinit
1134
1136
1137 if (space%is_periodic()) then
1138 assert(ubound(qq, 1) >= space%periodic_dim)
1139 assert(this%method == poisson_fft)
1140 end if
1141
1142 if (cam%omega > m_epsilon) then
1143 if (this%method /= poisson_fft) then
1144 write(message(1),'(a)') "Poisson solver with range separation is only implemented with FFT."
1145 call messages_fatal(1, namespace=namespace)
1146 end if
1147 end if
1148
1149 !We only reinitialize the poisson solver if needed
1150 reinit = .false.
1151 if (allocated(coulb%qq)) then
1152 reinit = any(abs(coulb%qq(1:space%periodic_dim) - qq(1:space%periodic_dim)) > m_epsilon)
1153 end if
1154 reinit = reinit .or. (abs(coulb%mu - cam%omega) > m_epsilon .and. cam%omega > m_epsilon)
1155 reinit = reinit .or. (abs(coulb%alpha - cam%alpha) > m_epsilon .and. cam%alpha > m_epsilon)
1156 reinit = reinit .or. (abs(coulb%beta - cam%beta) > m_epsilon .and. cam%beta > m_epsilon)
1157
1158 if (reinit) then
1159 !TODO: this should be a select case supporting other kernels.
1160 ! This means that we need an abstract object for kernels.
1161 select case (this%method)
1162 case (poisson_fft)
1163 ! Check that we are consistent: the Poisson solver supports must return 1 here
1164 assert(is_close(poisson_get_full_range_weight(this, cam), m_one))
1165
1166 call fourier_space_op_end(coulb)
1167
1168 safe_allocate(coulb%qq(1:space%dim))
1169 coulb%qq(1:space%periodic_dim) = qq(1:space%periodic_dim)
1170 coulb%qq(space%periodic_dim+1:space%dim) = vanishing_q
1171 !We must define the singularity if we specify a q vector and we do not use the short-range Coulomb potential
1172 coulb%singularity = optional_default(singul, m_zero)
1173 coulb%mu = cam%omega
1174 coulb%alpha = cam%alpha
1175 coulb%beta = cam%beta
1176 call poisson_fft_get_kernel(namespace, space, this%cube, coulb, this%kernel, &
1177 this%poisson_soft_coulomb_param)
1178 case default
1179 call messages_not_implemented("poisson_build_kernel with other methods than FFT", namespace=namespace)
1180 end select
1181 end if
1182
1183 pop_sub(poisson_build_kernel)
1184 end subroutine poisson_build_kernel
1185
1186 !----------------------------------------------------------------
1195 real(real64) function poisson_get_full_range_weight(this, cam) result(weight)
1196 type(poisson_t), intent(in) :: this
1197 type(xc_cam_t), intent(in) :: cam
1198
1199 select case (this%method)
1200 case (poisson_fft)
1201 weight = m_one
1202 case default
1203 if(cam%omega < m_epsilon) then
1204 weight = cam%alpha
1205 else if(cam%alpha > m_epsilon .and. cam%beta < m_epsilon) then
1206 weight = cam%alpha
1207 else if(cam%alpha < m_epsilon .and. cam%beta > m_epsilon) then
1208 weight = cam%beta
1209 else
1210 assert(.false.)
1211 end if
1212 end select
1214
1215#include "poisson_init_inc.F90"
1216#include "poisson_direct_inc.F90"
1217#include "poisson_direct_sm_inc.F90"
1218
1219#include "undef.F90"
1220#include "real.F90"
1221#include "poisson_inc.F90"
1222#include "undef.F90"
1223#include "complex.F90"
1224#include "poisson_inc.F90"
1225
1226end module poisson_oct_m
1227
1228!! Local Variables:
1229!! mode: f90
1230!! coding: utf-8
1231!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:181
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
This module implements batches of mesh functions.
Definition: batch.F90:135
This module handles the calculation mode.
integer, parameter, public p_strategy_kpoints
parallelization in k-points
subroutine, public cube_init(cube, nn, namespace, space, spacing, coord_system, fft_type, fft_library, dont_optimize, nn_out, mpi_grp, need_partition, tp_enlarge, blocksize)
Definition: cube.F90:204
subroutine, public cube_end(cube)
Definition: cube.F90:387
subroutine, public cube_init_cube_map(cube, mesh)
Definition: cube.F90:824
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:120
integer, parameter, public fft_none
global constants
Definition: fft.F90:174
integer, public fft_default_lib
Definition: fft.F90:250
integer, parameter, public fftlib_accel
Definition: fft.F90:179
integer, parameter, public fft_real
Definition: fft.F90:174
integer, parameter, public fft_complex
Definition: fft.F90:174
integer, parameter, public fftlib_pfft
Definition: fft.F90:179
integer, parameter, public fftlib_fftw
Definition: fft.F90:179
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the index, used for the mesh points.
Definition: index.F90:124
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
subroutine, public mesh_cube_parallel_map_end(this)
subroutine, public mesh_cube_parallel_map_init(this, mesh, cube)
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 mesh_double_box(space, mesh, alpha, db)
finds the dimension of a box doubled in the non-periodic dimensions
Definition: mesh.F90:285
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:897
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical pure function, public multicomm_strategy_is_parallel(mc, level)
Definition: multicomm.F90:719
logical pure function, public multicomm_have_slaves(this)
Definition: multicomm.F90:845
Some general things and nomenclature:
Definition: par_vec.F90:173
subroutine, public photon_mode_compute_dipoles(this, mesh)
Computes the polarization dipole.
subroutine, public photon_mode_add_poisson_terms(this, mesh, rho, pot)
subroutine, public photon_mode_end(this)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public photon_mode_init(this, namespace, dim, photon_free)
real(real64), public threshold
Definition: poisson_cg.F90:141
subroutine, public poisson_cg2(namespace, der, pot, rho)
Definition: poisson_cg.F90:217
subroutine, public poisson_cg1(namespace, der, corrector, pot, rho)
Definition: poisson_cg.F90:167
subroutine, public poisson_cg_init(thr, itr)
Definition: poisson_cg.F90:149
subroutine, public poisson_cg_end
Definition: poisson_cg.F90:161
subroutine, public poisson_corrections_end(this)
subroutine, public poisson_corrections_init(this, namespace, space, ml, mesh)
subroutine, public correct_rho(this, der, rho, rho_corrected, vh_correction)
integer, parameter, public poisson_fft_kernel_nocut
integer, parameter, public poisson_fft_kernel_cyl
subroutine, public zpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine, public poisson_fft_end(this)
subroutine, public poisson_fft_init(this, namespace, space, cube, kernel, soft_coulb_param, fullcube)
integer, parameter, public poisson_fft_kernel_pla
integer, parameter, public poisson_fft_kernel_none
integer, parameter, public poisson_fft_kernel_corrected
integer, parameter, public poisson_fft_kernel_sph
subroutine, public dpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine, public poisson_isf_end(this)
subroutine, public poisson_isf_init(this, namespace, mesh, cube, all_nodes_comm, init_world)
subroutine, public poisson_isf_solve(this, mesh, cube, pot, rho, all_nodes, sm)
subroutine, public poisson_multigrid_solver(this, namespace, der, pot, rho)
A multigrid Poisson solver with corrections at the boundaries.
subroutine, public poisson_multigrid_end(this)
subroutine, public poisson_no_solve(this, mesh, pot, rho)
Definition: poisson_no.F90:165
subroutine, public poisson_no_end(this)
Definition: poisson_no.F90:153
subroutine, public zpoisson_solve_sm(this, namespace, sm, pot, rho, all_nodes)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:2173
integer, parameter, public poisson_multigrid
Definition: poisson.F90:189
subroutine poisson_kernel_init(this, namespace, space, mc, stencil)
Definition: poisson.F90:1248
integer, parameter, public poisson_psolver
Definition: poisson.F90:189
subroutine, public dpoisson_solve_start(this, rho)
Definition: poisson.F90:1994
integer, parameter cmd_finish
Definition: poisson.F90:225
subroutine, public zpoisson_solve_finish(this, pot)
Definition: poisson.F90:2161
subroutine poisson_solve_direct(this, namespace, pot, rho)
Definition: poisson.F90:1452
integer, parameter, public poisson_fft
Definition: poisson.F90:189
subroutine, public zpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Definition: poisson.F90:791
subroutine, public poisson_solve_batch(this, namespace, potb, rhob, all_nodes, kernel)
Definition: poisson.F90:837
subroutine, public poisson_async_init(this, mc)
Definition: poisson.F90:1091
subroutine, public dpoisson_solve_sm(this, namespace, sm, pot, rho, all_nodes)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:2014
subroutine zpoisson_solve_real_and_imag_separately(this, namespace, pot, rho, all_nodes, kernel)
Definition: poisson.F90:742
subroutine, public poisson_init_sm(this, namespace, space, main, der, sm, method, force_cmplx)
Definition: poisson.F90:979
logical pure function poisson_solver_is_iterative(this)
Definition: poisson.F90:1084
subroutine, public poisson_slave_work(this, namespace)
Definition: poisson.F90:1119
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:872
integer, parameter cmd_poisson_solve
Definition: poisson.F90:225
integer, parameter, public poisson_cg
Definition: poisson.F90:189
subroutine, public poisson_build_kernel(this, namespace, space, coulb, qq, cam, singul)
Definition: poisson.F90:1136
subroutine, public dpoisson_solve_finish(this, pot)
Definition: poisson.F90:2002
subroutine, public poisson_init(this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx)
Definition: poisson.F90:233
subroutine, public zpoisson_solve_start(this, rho)
Definition: poisson.F90:2153
subroutine, public poisson_async_end(this, mc)
Definition: poisson.F90:1103
integer, parameter, public poisson_cg_corrected
Definition: poisson.F90:189
integer, parameter, public poisson_isf
Definition: poisson.F90:189
real(real64) function, public poisson_get_full_range_weight(this, cam)
Most Poisson solvers do not implement Coulomb attenuated potentials, and can only be used for global ...
Definition: poisson.F90:1208
integer, parameter, public poisson_null
Definition: poisson.F90:189
integer, parameter, public poisson_no
Definition: poisson.F90:189
logical pure function, public poisson_is_async(this)
Definition: poisson.F90:1127
subroutine, public poisson_end(this)
Definition: poisson.F90:689
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_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
This module defines stencils used in Octopus.
Definition: stencil.F90:137
subroutine, public submesh_init_cube_map(sm, space)
Definition: submesh.F90:932
subroutine, public submesh_get_cube_dim(sm, space, db)
finds the dimension of a box containing the submesh
Definition: submesh.F90:897
type(type_t), public type_float
Definition: types.F90:135
This module defines the unit system, used for input and output.
Class defining batches of mesh functions.
Definition: batch.F90:161
Class implementing a box that is a union of spheres. We do this in a specific class instead of using ...
class representing derivatives
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)