Octopus
multigrid.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade
2!! Copyright (C) 2021 S. Ohlmann
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 multigrid_oct_m
23 use batch_oct_m
26 use debug_oct_m
28 use global_oct_m
29 use index_oct_m
31 use math_oct_m
32 use mesh_oct_m
35 use mpi_oct_m
38 use parser_oct_m
40 use space_oct_m
44
45 implicit none
46
47 private
48 public :: &
65
66 integer, parameter, public :: &
67 INJECTION = 1, &
68 fullweight = 2
69
71 ! Components are public by default
72 type(transfer_table_t) :: tt
73 type(mesh_t), pointer :: mesh => null()
74 type(derivatives_t), pointer :: der => null()
75 end type multigrid_level_t
76
77 type multigrid_t
78 private
79 integer, public :: n_levels
80 type(multigrid_level_t), allocatable, public :: level(:)
81
82 integer :: tp
83 integer, allocatable :: sp(:)
84 integer, allocatable :: ep(:)
85 integer, allocatable :: ep_part(:)
86 end type multigrid_t
87
88
89contains
90
91 ! ---------------------------------------------------------
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
101
102 integer :: i, n_levels, np, order
103
104 push_sub(multigrid_init)
105
106 if (.not. present(nlevels)) then
107 !%Variable MultigridLevels
108 !%Type integer
109 !%Default 4
110 !%Section Mesh
111 !%Description
112 !% Number of levels in the grid hierarchy used for multigrid. Positive
113 !% numbers indicate an absolute number of levels, negative
114 !% numbers are subtracted from the maximum number of levels possible.
115 !%Option max_levels 0
116 !% Calculate the optimal number of levels for the grid.
117 !%End
118
119 call parse_variable(namespace, 'MultigridLevels', 4, n_levels)
120
121 ! default:
122 order = der%order
123 else
124 n_levels = nlevels
125 write(message(1), '(a,i1,a)') "Set number of multigrid levels to ", n_levels, ". This ignores the value of MultigridLevels."
126 call messages_info(1, namespace=namespace)
127
128 !%Variable MultigridDerivativesOrder
129 !%Type integer
130 !%Default 1
131 !%Section Mesh::Derivatives
132 !%Description
133 !% This variable gives the discretization order for the approximation of
134 !% the differential operators on the different levels of the multigrid.
135 !% For more details, see the variable DerivativesOrder.
136 !% For star-general stencils, the minimum is set to 2.
137 !%End
138 call parse_variable(namespace, 'MultigridDerivativesOrder', 1, order)
139 ! set order to a minimum of 2 for general star stencil, fails otherwise
140 if (der%stencil_type == der_stargeneral) then
141 order = max(2, order)
142 end if
143 end if
144
145 ! require at least one cell per direction on the coarsest grid and adapt n_levels if needed
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."
148 call messages_warning(1)
149 do i = n_levels, 1, -1
150 if (all(mesh%idx%ll >= 2**i)) then
151 n_levels = i
152 exit
153 end if
154 end do
155 end if
156
157 if (n_levels <= 0) then
158 n_levels = n_levels - 3
159 np = mesh%np
160
161 do while (np > 1)
162 np = np / 8
163 n_levels = n_levels + 1
164 end do
165 else
166 n_levels = n_levels - 1
167 end if
169 mgrid%n_levels = n_levels
170
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
180 call messages_info(1, namespace=namespace)
181
182 do i = 1, mgrid%n_levels
183 safe_allocate(mgrid%level(i)%mesh)
184 safe_allocate(mgrid%level(i)%der)
185
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)
189
190 call mesh_init_stage_3(mgrid%level(i)%mesh, namespace, space, stencil, mc, parent = mgrid%level(i - 1)%mesh)
191
192 call multigrid_get_transfer_tables(mgrid%level(i)%tt, space, mgrid%level(i-1)%mesh, mgrid%level(i)%mesh)
193
194 call derivatives_build(mgrid%level(i)%der, namespace, space, mgrid%level(i)%mesh)
195
196 call mesh_write_info(mgrid%level(i)%mesh, namespace=namespace)
197
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
202 end do
203
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))
207
208 mgrid%tp = 0
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
214 end do
215
216 pop_sub(multigrid_init)
217 end subroutine multigrid_init
218
219 ! ---------------------------------------------------------
221 subroutine multigrid_get_transfer_tables(tt, space, fine, coarse)
222 type(transfer_table_t), intent(inout) :: tt
223 class(space_t), intent(in) :: space
224 type(mesh_t), intent(in) :: fine, coarse
225
226 integer :: i, i1, i2, i4, i8, pt
227 integer :: x(space%dim), mod2(space%dim), idx(space%dim)
228
230
231 ! Currently this routine only works for 3D
232 assert(space%dim == 3)
233
234 tt%n_coarse = coarse%np
235 safe_allocate(tt%to_coarse(1:tt%n_coarse))
236
237 ! GENERATE THE TABLE TO MAP FROM THE FINE TO THE COARSE GRID
238 do i = 1, tt%n_coarse
239 ! locate the equivalent fine grid point
240 call mesh_local_index_to_coords(coarse, i, idx)
241 tt%to_coarse(i) = mesh_local_index_from_coords(fine, 2*idx)
242 end do
243
244 ! count
245 tt%n_fine = fine%np
246 safe_allocate(tt%fine_i(1:tt%n_fine))
247
248 tt%n_fine1 = 0
249 tt%n_fine2 = 0
250 tt%n_fine4 = 0
251 tt%n_fine8 = 0
252 do i = 1, tt%n_fine
253 call mesh_local_index_to_coords(fine, i, idx)
254 mod2 = mod(idx, 2)
255
256 pt = sum(abs(mod2))
257
258 select case (pt)
259 case (0)
260 tt%n_fine1 = tt%n_fine1 + 1
261 tt%fine_i(i) = 1
262 case (1)
263 tt%n_fine2 = tt%n_fine2 + 1
264 tt%fine_i(i) = 2
265 case (2)
266 tt%n_fine4 = tt%n_fine4 + 1
267 tt%fine_i(i) = 4
268 case (3)
269 tt%n_fine8 = tt%n_fine8 + 1
270 tt%fine_i(i) = 8
271 end select
272 end do
273
274 assert(tt%n_fine1 + tt%n_fine2 + tt%n_fine4 + tt%n_fine8 == tt%n_fine)
275
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))
280
281 ! and now build the tables
282 i1 = 0
283 i2 = 0
284 i4 = 0
285 i8 = 0
286 do i = 1, fine%np
287 call mesh_local_index_to_coords(fine, i, idx)
288 x = idx/2
289 mod2 = mod(idx, 2)
290
291 pt = sum(abs(mod2))
292
293 select case (pt)
294 case (0)
295 i1 = i1 + 1
296 tt%to_fine1(1, i1) = mesh_local_index_from_coords(coarse, x)
297
298 case (1)
299 i2 = i2 + 1
300 tt%to_fine2(1, i2) = mesh_local_index_from_coords(coarse, x)
301 tt%to_fine2(2, i2) = mesh_local_index_from_coords(coarse, x + mod2)
302
303 case (2)
304 i4 = i4 + 1
305 tt%to_fine4(1, i4) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) + mod2(3)])
306 tt%to_fine4(2, i4) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) + mod2(3)])
307 tt%to_fine4(3, i4) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) + mod2(2), x(3) ])
308 tt%to_fine4(4, i4) = mesh_local_index_from_coords(coarse, x + mod2)
309
310 case (3)
311 i8 = i8 + 1
312 tt%to_fine8(1, i8) = mesh_local_index_from_coords(coarse, x)
313 tt%to_fine8(2, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) ])
314 tt%to_fine8(3, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) ])
315 tt%to_fine8(4, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) , x(3) + mod2(3)])
316 tt%to_fine8(5, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) + mod2(3)])
317 tt%to_fine8(6, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) + mod2(3)])
318 tt%to_fine8(7, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) + mod2(2), x(3) ])
319 tt%to_fine8(8, i8) = mesh_local_index_from_coords(coarse, x + mod2)
320
321 end select
322
323 end do
324
325 assert(i1 == tt%n_fine1 .and. i2 == tt%n_fine2 .and. i4 == tt%n_fine4 .and. i8 == tt%n_fine8)
326
327
329 end subroutine multigrid_get_transfer_tables
330
331 !---------------------------------------------------------------------------------
334 !---------------------------------------------------------------------------------
335 subroutine multigrid_mesh_half(space, namespace, mesh_in, mesh_out, stencil, grp)
336 class(space_t), intent(in) :: space
337 type(namespace_t), intent(in) :: namespace
338 type(mesh_t), target, intent(in) :: mesh_in
339 type(mesh_t), intent(inout) :: mesh_out
340 type(stencil_t), intent(in) :: stencil
341 type(mpi_grp_t), intent(in) :: grp
342
343 integer :: idim
344
345 push_sub(multigrid_mesh_half)
346
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
352
353 safe_allocate(mesh_out%spacing(1:space%dim))
354 mesh_out%spacing(:) = 2*mesh_in%spacing(:)
355
356 call index_init(mesh_out%idx, space%dim)
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
361
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))
366 end do
367
368 call mesh_init_stage_2(mesh_out, namespace, space, mesh_out%box, stencil, grp)
369
370 pop_sub(multigrid_mesh_half)
371 end subroutine multigrid_mesh_half
372
373 !---------------------------------------------------------------------------------
374 subroutine multigrid_mesh_double(space, namespace, mesh_in, mesh_out, stencil, grp)
375 class(space_t), intent(in) :: space
376 type(namespace_t), intent(in) :: namespace
377 type(mesh_t), target, intent(in) :: mesh_in
378 type(mesh_t), intent(inout) :: mesh_out
379 type(stencil_t), intent(in) :: stencil
380 type(mpi_grp_t), intent(in) :: grp
381
382 integer :: idim
383 push_sub(multigrid_mesh_double)
384
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
390
391 safe_allocate(mesh_out%spacing(1:space%dim))
392 mesh_out%spacing(:) = m_half*mesh_in%spacing(:)
393
394 call index_init(mesh_out%idx, space%dim)
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
398 ! We need to make the possible number of grid points larger by one for
399 ! the non-periodic dimensions because the spacing is only half of the
400 ! original mesh and thus we could get points in the new boundary
401 ! that are still inside the simulation box.
402 ! For the periodic dimensions, we are anyway commensurate with the size
403 ! of the box, so we are still commensurate when taking twice the number
404 ! of points.
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
408 end do
409 mesh_out%idx%ll(:) = mesh_out%idx%nr(2, :) - mesh_out%idx%nr(1, :) + 1
410
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))
415 end do
416
417 call mesh_init_stage_2(mesh_out, namespace, space, mesh_out%box, stencil, grp)
418
419 pop_sub(multigrid_mesh_double)
420 end subroutine multigrid_mesh_double
421
422 ! ---------------------------------------------------------
423 subroutine multigrid_end(mgrid)
424 type(multigrid_t), target, intent(inout) :: mgrid
425
426 integer :: i
427 type(multigrid_level_t), pointer :: level
428
429 push_sub(multigrid_end)
431 safe_deallocate_a(mgrid%sp)
432 safe_deallocate_a(mgrid%ep)
433 safe_deallocate_a(mgrid%ep_part)
434
435 safe_deallocate_a(mgrid%level(0)%tt%fine_i)
436
437 do i = 1, mgrid%n_levels
438 level => mgrid%level(i)
439
440 call derivatives_end(level%der)
441 call mesh_end(level%mesh)
442 safe_deallocate_p(level%mesh)
443 safe_deallocate_p(level%der)
444
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)
451 end do
452
453 safe_deallocate_a(mgrid%level)
454
455 pop_sub(multigrid_end)
456 end subroutine multigrid_end
457
458 !---------------------------------------------------------------------------------
459 integer function multigrid_number_of_levels(base_der) result(number)
460 type(derivatives_t), target, intent(in) :: base_der
461
462 type(derivatives_t), pointer :: next_der
463
464 next_der => base_der%coarser
465
466 number = 0
467 next_der => base_der%coarser
468 do while (associated(next_der))
469 number = number + 1
470 next_der => next_der%coarser
471 end do
472
473 end function multigrid_number_of_levels
474
475 !---------------------------------------------------------------------------------
476 subroutine multigrid_build_stencil(dim, weight, shift)
477 integer, intent(in) :: dim
478 real(real64), intent(inout) :: weight(:)
479 integer, intent(inout) :: shift(:,:)
480
481 integer :: nn, di, dj, dk, dd
482
484
485 assert(ubound(weight, dim=1) == 3**dim)
486 assert(ubound(shift, dim=1) == dim)
487 assert(ubound(shift, dim=2) == 3**dim)
488
489 nn = 0
490 select case(dim)
491 case(1)
492 do di = -1, 1
493 dd = abs(di)
494 nn = nn + 1
495 weight(nn) = m_half**dd
496 shift(1,nn) = di
497 end do
498 case(2)
499 do di = -1, 1
500 do dj = -1, 1
501 dd = abs(di)+abs(dj)
502 nn = nn + 1
503 weight(nn) = m_half**dd
504 shift(1,nn) = di
505 shift(2,nn) = dj
506 end do
507 end do
508 case(3)
509 do di = -1, 1
510 do dj = -1, 1
511 do dk = -1, 1
512 dd = abs(di)+abs(dj)+abs(dk)
513 nn = nn + 1
514 weight(nn) = m_half**dd
515 shift(1,nn) = di
516 shift(2,nn) = dj
517 shift(3,nn) = dk
518 end do
519 end do
520 end do
521 end select
522
524 end subroutine multigrid_build_stencil
525
526
527#include "undef.F90"
528#include "real.F90"
529#include "multigrid_inc.F90"
530
531#include "undef.F90"
532#include "complex.F90"
533#include "multigrid_inc.F90"
534
535end module multigrid_oct_m
536
537!! Local Variables:
538!! mode: f90
539!! coding: utf-8
540!! End:
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
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
Definition: global.F90:206
This module implements the index, used for the mesh points.
Definition: index.F90:124
subroutine, public index_init(idx, dim)
This subroutine allocates memory and initializes some components.
Definition: index.F90:220
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module contains subroutines, related to the initialization of the mesh.
Definition: mesh_init.F90:119
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....
Definition: mesh_init.F90:531
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.
Definition: mesh_init.F90:292
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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.
Definition: mesh.F90:938
subroutine, public mesh_write_info(this, iunit, namespace)
Definition: mesh.F90:310
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.
Definition: mesh.F90:950
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:686
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public dmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
Definition: multigrid.F90:989
subroutine, public dmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, set_bc)
Definition: multigrid.F90:691
integer function, public multigrid_number_of_levels(base_der)
Definition: multigrid.F90:555
subroutine, public multigrid_end(mgrid)
Definition: multigrid.F90:519
subroutine, public dmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
Definition: multigrid.F90:770
subroutine, public multigrid_mesh_double(space, namespace, mesh_in, mesh_out, stencil, grp)
Definition: multigrid.F90:470
integer, parameter, public fullweight
Definition: multigrid.F90:161
subroutine, public zmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
Definition: multigrid.F90:1499
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...
Definition: multigrid.F90:431
subroutine, public zmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, set_bc)
Definition: multigrid.F90:1201
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
Definition: multigrid.F90:188
subroutine, public zmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
Definition: multigrid.F90:1280
subroutine, public dmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb)
Definition: multigrid.F90:885
subroutine multigrid_build_stencil(dim, weight, shift)
Definition: multigrid.F90:572
subroutine, public multigrid_get_transfer_tables(tt, space, fine, coarse)
creates the lookup tables to go between the coarse and fine meshes
Definition: multigrid.F90:317
subroutine, public zmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb)
Definition: multigrid.F90:1395
Some general things and nomenclature:
Definition: par_vec.F90:173
This module defines stencils used in Octopus.
Definition: stencil.F90:137
class representing derivatives
Describes mesh distribution to nodes.
Definition: mesh.F90:187
This is defined even when running serial.
Definition: mpi.F90:144
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:165