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