66 integer,
parameter,
public :: &
72 type(transfer_table_t) :: tt
73 type(mesh_t),
pointer :: mesh => null()
74 type(derivatives_t),
pointer :: der => null()
79 integer,
public :: n_levels
80 type(multigrid_level_t),
allocatable,
public :: level(:)
83 integer,
allocatable :: sp(:)
84 integer,
allocatable :: ep(:)
85 integer,
allocatable :: ep_part(:)
92 subroutine multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
93 type(multigrid_t),
target,
intent(out) :: mgrid
94 type(namespace_t),
intent(in) :: namespace
95 class(space_t),
intent(in) :: space
96 class(mesh_t),
target,
intent(in) :: mesh
97 type(derivatives_t),
target,
intent(in) :: der
98 type(stencil_t),
intent(in) :: stencil
99 type(multicomm_t),
intent(in) :: mc
100 integer,
optional,
intent(in) :: nlevels
102 integer :: i, n_levels, np, order
106 if (.not.
present(nlevels))
then
125 write(
message(1),
'(a,i1,a)')
"Set number of multigrid levels to ", n_levels,
". This ignores the value of MultigridLevels."
138 call parse_variable(namespace,
'MultigridDerivativesOrder', 1, order)
141 order = max(2, order)
146 if (any(mesh%idx%ll < 2**n_levels))
then
147 message(1) =
"Grid too small for the multigrid solver. Reducing the number of levels."
149 do i = n_levels, 1, -1
150 if (all(mesh%idx%ll >= 2**i))
then
157 if (n_levels <= 0)
then
158 n_levels = n_levels - 3
163 n_levels = n_levels + 1
166 n_levels = n_levels - 1
169 mgrid%n_levels = n_levels
171 safe_allocate(mgrid%level(0:n_levels))
173 mgrid%level(0)%mesh => mesh
174 mgrid%level(0)%der => der
176 mgrid%level(0)%tt%n_fine = mesh%np
177 safe_allocate(mgrid%level(0)%tt%fine_i(1:mesh%np))
179 write(
message(1),
'(a,i3)')
"Multigrid levels:", n_levels + 1
182 do i = 1, mgrid%n_levels
183 safe_allocate(mgrid%level(i)%mesh)
184 safe_allocate(mgrid%level(i)%der)
186 call multigrid_mesh_half(space, namespace, mgrid%level(i-1)%mesh, mgrid%level(i)%mesh, stencil, mc%base_grp)
188 call derivatives_init(mgrid%level(i)%der, namespace, space, mesh%coord_system, order=order)
190 call mesh_init_stage_3(mgrid%level(i)%mesh, namespace, space, stencil, mc, parent = mgrid%level(i - 1)%mesh)
194 call derivatives_build(mgrid%level(i)%der, namespace, space, mgrid%level(i)%mesh)
198 mgrid%level(i)%der%finer => mgrid%level(i - 1)%der
199 mgrid%level(i - 1)%der%coarser => mgrid%level(i)%der
200 mgrid%level(i)%der%to_finer => mgrid%level(i)%tt
201 mgrid%level(i - 1)%der%to_coarser => mgrid%level(i)%tt
204 safe_allocate(mgrid%sp(0:mgrid%n_levels))
205 safe_allocate(mgrid%ep(0:mgrid%n_levels))
206 safe_allocate(mgrid%ep_part(0:mgrid%n_levels))
209 do i = 0, mgrid%n_levels
210 mgrid%sp(i) = 1 + mgrid%tp
211 mgrid%ep(i) = mgrid%tp + mgrid%level(i)%mesh%np
212 mgrid%tp = mgrid%tp + mgrid%level(i)%mesh%np_part
213 mgrid%ep_part(i) = mgrid%tp
223 class(
space_t),
intent(in) :: space
224 type(
mesh_t),
intent(in) :: fine, coarse
226 integer :: i, i1, i2, i4, i8, pt
227 integer :: x(space%dim), mod2(space%dim), idx(space%dim)
232 assert(space%dim == 3)
234 tt%n_coarse = coarse%np
235 safe_allocate(tt%to_coarse(1:tt%n_coarse))
238 do i = 1, tt%n_coarse
246 safe_allocate(tt%fine_i(1:tt%n_fine))
260 tt%n_fine1 = tt%n_fine1 + 1
263 tt%n_fine2 = tt%n_fine2 + 1
266 tt%n_fine4 = tt%n_fine4 + 1
269 tt%n_fine8 = tt%n_fine8 + 1
274 assert(tt%n_fine1 + tt%n_fine2 + tt%n_fine4 + tt%n_fine8 == tt%n_fine)
276 safe_allocate(tt%to_fine1(1, 1:tt%n_fine1))
277 safe_allocate(tt%to_fine2(1:2, 1:tt%n_fine2))
278 safe_allocate(tt%to_fine4(1:4, 1:tt%n_fine4))
279 safe_allocate(tt%to_fine8(1:8, 1:tt%n_fine8))
325 assert(i1 == tt%n_fine1 .and. i2 == tt%n_fine2 .and. i4 == tt%n_fine4 .and. i8 == tt%n_fine8)
336 class(
space_t),
intent(in) :: space
338 type(
mesh_t),
target,
intent(in) :: mesh_in
339 type(
mesh_t),
intent(inout) :: mesh_out
347 mesh_out%box => mesh_in%box
348 mesh_out%idx%dim = mesh_in%idx%dim
349 mesh_out%use_curvilinear = mesh_in%use_curvilinear
350 mesh_out%masked_periodic_boundaries = mesh_in%masked_periodic_boundaries
351 mesh_out%coord_system => mesh_in%coord_system
353 safe_allocate(mesh_out%spacing(1:space%dim))
354 mesh_out%spacing(:) = 2*mesh_in%spacing(:)
357 mesh_out%idx%enlarge(:) = mesh_in%idx%enlarge(:)
358 mesh_out%idx%nr(1,:) = (mesh_in%idx%nr(1,:)+mesh_in%idx%enlarge(:))/2
359 mesh_out%idx%nr(2,:) = (mesh_in%idx%nr(2,:)-mesh_in%idx%enlarge(:))/2
360 mesh_out%idx%ll(:) = mesh_out%idx%nr(2, :) - mesh_out%idx%nr(1, :) + 1
362 mesh_out%idx%stride(1) = 1
363 do idim = 2, mesh_out%idx%dim + 1
364 mesh_out%idx%stride(idim) = mesh_out%idx%stride(idim-1) * &
365 (mesh_out%idx%ll(idim-1) + 2*mesh_out%idx%enlarge(idim-1))
375 class(
space_t),
intent(in) :: space
377 type(
mesh_t),
target,
intent(in) :: mesh_in
378 type(
mesh_t),
intent(inout) :: mesh_out
385 mesh_out%box => mesh_in%box
386 mesh_out%idx%dim = mesh_in%idx%dim
387 mesh_out%use_curvilinear = mesh_in%use_curvilinear
388 mesh_out%masked_periodic_boundaries = mesh_in%masked_periodic_boundaries
389 mesh_out%coord_system => mesh_in%coord_system
391 safe_allocate(mesh_out%spacing(1:space%dim))
392 mesh_out%spacing(:) =
m_half*mesh_in%spacing(:)
395 mesh_out%idx%enlarge(:) = mesh_in%idx%enlarge(:)
396 mesh_out%idx%nr(1,:) = (mesh_in%idx%nr(1,:)+mesh_in%idx%enlarge(:))*2
397 mesh_out%idx%nr(2,:) = (mesh_in%idx%nr(2,:)-mesh_in%idx%enlarge(:))*2
405 do idim = space%periodic_dim + 1, space%dim
406 mesh_out%idx%nr(1, idim) = mesh_out%idx%nr(1, idim) - 1
407 mesh_out%idx%nr(2, idim) = mesh_out%idx%nr(2, idim) + 1
409 mesh_out%idx%ll(:) = mesh_out%idx%nr(2, :) - mesh_out%idx%nr(1, :) + 1
411 mesh_out%idx%stride(1) = 1
412 do idim = 2, space%dim+1
413 mesh_out%idx%stride(idim) = mesh_out%idx%stride(idim-1) * &
414 (mesh_out%idx%ll(idim-1) + 2*mesh_out%idx%enlarge(idim-1))
431 safe_deallocate_a(mgrid%sp)
432 safe_deallocate_a(mgrid%ep)
433 safe_deallocate_a(mgrid%ep_part)
435 safe_deallocate_a(mgrid%level(0)%tt%fine_i)
437 do i = 1, mgrid%n_levels
438 level => mgrid%level(i)
442 safe_deallocate_p(level%mesh)
443 safe_deallocate_p(level%der)
445 safe_deallocate_a(level%tt%to_coarse)
446 safe_deallocate_a(level%tt%to_fine1)
447 safe_deallocate_a(level%tt%to_fine2)
448 safe_deallocate_a(level%tt%to_fine4)
449 safe_deallocate_a(level%tt%to_fine8)
450 safe_deallocate_a(level%tt%fine_i)
453 safe_deallocate_a(mgrid%level)
464 next_der => base_der%coarser
467 next_der => base_der%coarser
468 do while (
associated(next_der))
470 next_der => next_der%coarser
477 integer,
intent(in) :: dim
478 real(real64),
intent(inout) :: weight(:)
479 integer,
intent(inout) :: shift(:,:)
481 integer :: nn, di, dj, dk, dd
485 assert(ubound(weight, dim=1) == 3**dim)
486 assert(ubound(shift, dim=1) == dim)
487 assert(ubound(shift, dim=2) == 3**dim)
512 dd = abs(di)+abs(dj)+abs(dk)
529#include "multigrid_inc.F90"
532#include "complex.F90"
533#include "multigrid_inc.F90"
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
build the derivatives object:
subroutine, public derivatives_init(der, namespace, space, coord_system, order)
subroutine, public derivatives_end(der)
integer, parameter, public der_stargeneral
real(real64), parameter, public m_half
This module implements the index, used for the mesh points.
subroutine, public index_init(idx, dim)
This subroutine allocates memory and initializes some components.
This module is intended to contain "only mathematical" functions and procedures.
This module contains subroutines, related to the initialization of the mesh.
subroutine, public mesh_init_stage_3(mesh, namespace, space, stencil, mc, parent, regenerate)
When running parallel in domains, stencil and np_stencil are needed to compute the ghost points....
subroutine, public mesh_init_stage_2(mesh, namespace, space, box, stencil, grp, regenerate)
This subroutine creates the global array of spatial indices and the inverse mapping.
This module defines the meshes, which are used in Octopus.
integer function, public mesh_local_index_from_coords(mesh, ix)
This function returns the local index of the point for a given vector of integer coordinates.
subroutine, public mesh_write_info(this, iunit, namespace)
subroutine, public mesh_local_index_to_coords(mesh, ip, ix)
Given a local point index, this function returns the set of integer coordinates of the point.
recursive subroutine, public mesh_end(this)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module handles the communicators for the various parallelization strategies.
subroutine, public dmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
subroutine, public dmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, set_bc)
integer function, public multigrid_number_of_levels(base_der)
subroutine, public multigrid_end(mgrid)
subroutine, public dmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
subroutine, public multigrid_mesh_double(space, namespace, mesh_in, mesh_out, stencil, grp)
integer, parameter, public fullweight
subroutine, public zmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
subroutine, public multigrid_mesh_half(space, namespace, mesh_in, mesh_out, stencil, grp)
Creates a mesh that has twice the spacing betwen the points than the in mesh. This is used in the mul...
subroutine, public zmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, set_bc)
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
subroutine, public zmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
subroutine, public dmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb)
subroutine multigrid_build_stencil(dim, weight, shift)
subroutine, public multigrid_get_transfer_tables(tt, space, fine, coarse)
creates the lookup tables to go between the coarse and fine meshes
subroutine, public zmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb)
Some general things and nomenclature:
This module defines stencils used in Octopus.
class representing derivatives
Describes mesh distribution to nodes.
This is defined even when running serial.
The class representing the stencil, which is used for non-local mesh operations.