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
30 use math_oct_m
31 use mesh_oct_m
36 use parser_oct_m
38 use space_oct_m
42
43 implicit none
44
45 private
46 public :: &
63
64 integer, parameter, public :: &
65 INJECTION = 1, &
66 fullweight = 2
67
69 ! Components are public by default
70 type(transfer_table_t) :: tt
71 type(mesh_t), pointer :: mesh => null()
72 type(derivatives_t), pointer :: der => null()
73 end type multigrid_level_t
74
75 type multigrid_t
76 private
77 integer, public :: n_levels
78 type(multigrid_level_t), allocatable, public :: level(:)
79
80 integer :: tp
81 integer, allocatable :: sp(:)
82 integer, allocatable :: ep(:)
83 integer, allocatable :: ep_part(:)
84 end type multigrid_t
85
86
87contains
88
89 ! ---------------------------------------------------------
90 subroutine multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
91 type(multigrid_t), target, intent(out) :: mgrid
92 type(namespace_t), intent(in) :: namespace
93 class(space_t), intent(in) :: space
94 class(mesh_t), target, intent(in) :: mesh
95 type(derivatives_t), target, intent(in) :: der
96 type(stencil_t), intent(in) :: stencil
97 type(multicomm_t), intent(in) :: mc
98 integer, optional, intent(in) :: nlevels
99
100 integer :: i, n_levels, np, order
101
102 push_sub(multigrid_init)
103
104 if (.not. present(nlevels)) then
105 !%Variable MultigridLevels
106 !%Type integer
107 !%Default 4
108 !%Section Mesh
109 !%Description
110 !% Number of levels in the grid hierarchy used for multigrid. Positive
111 !% numbers indicate an absolute number of levels, negative
112 !% numbers are subtracted from the maximum number of levels possible.
113 !%Option max_levels 0
114 !% Calculate the optimal number of levels for the grid.
115 !%End
116
117 call parse_variable(namespace, 'MultigridLevels', 4, n_levels)
118
119 ! default:
120 order = der%order
121 else
122 n_levels = nlevels
123 write(message(1), '(a,i1,a)') "Set number of multigrid levels to ", n_levels, ". This ignores the value of MultigridLevels."
124 call messages_info(1, namespace=namespace)
125
126 !%Variable MultigridDerivativesOrder
127 !%Type integer
128 !%Default 1
129 !%Section Mesh::Derivatives
130 !%Description
131 !% This variable gives the discretization order for the approximation of
132 !% the differential operators on the different levels of the multigrid.
133 !% For more details, see the variable DerivativesOrder.
134 !% For star-general stencils, the minimum is set to 2.
135 !%End
136 call parse_variable(namespace, 'MultigridDerivativesOrder', 1, order)
137 ! set order to a minimum of 2 for general star stencil, fails otherwise
138 if (der%stencil_type == der_stargeneral) then
139 order = max(2, order)
140 end if
141 end if
142
143 if (n_levels <= 0) then
144 n_levels = n_levels - 3
145 np = mesh%np
146
147 do while (np > 1)
148 np = np / 8
149 n_levels = n_levels + 1
150 end do
151 else
152 n_levels = n_levels - 1
153 end if
154
155 mgrid%n_levels = n_levels
156
157 safe_allocate(mgrid%level(0:n_levels))
158
159 mgrid%level(0)%mesh => mesh
160 mgrid%level(0)%der => der
162 mgrid%level(0)%tt%n_fine = mesh%np
163 safe_allocate(mgrid%level(0)%tt%fine_i(1:mesh%np))
165 write(message(1), '(a,i3)') "Multigrid levels:", n_levels + 1
166 call messages_info(1, namespace=namespace)
167
168 do i = 1, mgrid%n_levels
169 safe_allocate(mgrid%level(i)%mesh)
170 safe_allocate(mgrid%level(i)%der)
172 call multigrid_mesh_half(space, namespace, mgrid%level(i-1)%mesh, mgrid%level(i)%mesh, stencil)
174 call derivatives_init(mgrid%level(i)%der, namespace, space, mesh%coord_system, order=order)
176 call mesh_init_stage_3(mgrid%level(i)%mesh, namespace, space, stencil, mc, parent = mgrid%level(i - 1)%mesh)
177
178 call multigrid_get_transfer_tables(mgrid%level(i)%tt, space, mgrid%level(i-1)%mesh, mgrid%level(i)%mesh)
179
180 call derivatives_build(mgrid%level(i)%der, namespace, space, mgrid%level(i)%mesh)
181
182 call mesh_write_info(mgrid%level(i)%mesh, namespace=namespace)
184 mgrid%level(i)%der%finer => mgrid%level(i - 1)%der
185 mgrid%level(i - 1)%der%coarser => mgrid%level(i)%der
186 mgrid%level(i)%der%to_finer => mgrid%level(i)%tt
187 mgrid%level(i - 1)%der%to_coarser => mgrid%level(i)%tt
188 end do
189
190 safe_allocate(mgrid%sp(0:mgrid%n_levels))
191 safe_allocate(mgrid%ep(0:mgrid%n_levels))
192 safe_allocate(mgrid%ep_part(0:mgrid%n_levels))
193
194 mgrid%tp = 0
195 do i = 0, mgrid%n_levels
196 mgrid%sp(i) = 1 + mgrid%tp
197 mgrid%ep(i) = mgrid%tp + mgrid%level(i)%mesh%np
198 mgrid%tp = mgrid%tp + mgrid%level(i)%mesh%np_part
199 mgrid%ep_part(i) = mgrid%tp
200 end do
201
202 pop_sub(multigrid_init)
203 end subroutine multigrid_init
204
205 ! ---------------------------------------------------------
207 subroutine multigrid_get_transfer_tables(tt, space, fine, coarse)
208 type(transfer_table_t), intent(inout) :: tt
209 class(space_t), intent(in) :: space
210 type(mesh_t), intent(in) :: fine, coarse
211
212 integer :: i, i1, i2, i4, i8, pt
213 integer :: x(space%dim), mod2(space%dim), idx(space%dim)
214
216
217 ! Currently this routine only works for 3D
218 assert(space%dim == 3)
219
220 tt%n_coarse = coarse%np
221 safe_allocate(tt%to_coarse(1:tt%n_coarse))
222
223 ! GENERATE THE TABLE TO MAP FROM THE FINE TO THE COARSE GRID
224 do i = 1, tt%n_coarse
225 ! locate the equivalent fine grid point
226 call mesh_local_index_to_coords(coarse, i, idx)
227 tt%to_coarse(i) = mesh_local_index_from_coords(fine, 2*idx)
228 end do
229
230 ! count
231 tt%n_fine = fine%np
232 safe_allocate(tt%fine_i(1:tt%n_fine))
233
234 tt%n_fine1 = 0
235 tt%n_fine2 = 0
236 tt%n_fine4 = 0
237 tt%n_fine8 = 0
238 do i = 1, tt%n_fine
239 call mesh_local_index_to_coords(fine, i, idx)
240 mod2 = mod(idx, 2)
241
242 pt = sum(abs(mod2))
243
244 select case (pt)
245 case (0)
246 tt%n_fine1 = tt%n_fine1 + 1
247 tt%fine_i(i) = 1
248 case (1)
249 tt%n_fine2 = tt%n_fine2 + 1
250 tt%fine_i(i) = 2
251 case (2)
252 tt%n_fine4 = tt%n_fine4 + 1
253 tt%fine_i(i) = 4
254 case (3)
255 tt%n_fine8 = tt%n_fine8 + 1
256 tt%fine_i(i) = 8
257 end select
258 end do
259
260 assert(tt%n_fine1 + tt%n_fine2 + tt%n_fine4 + tt%n_fine8 == tt%n_fine)
261
262 safe_allocate(tt%to_fine1(1, 1:tt%n_fine1))
263 safe_allocate(tt%to_fine2(1:2, 1:tt%n_fine2))
264 safe_allocate(tt%to_fine4(1:4, 1:tt%n_fine4))
265 safe_allocate(tt%to_fine8(1:8, 1:tt%n_fine8))
266
267 ! and now build the tables
268 i1 = 0
269 i2 = 0
270 i4 = 0
271 i8 = 0
272 do i = 1, fine%np
273 call mesh_local_index_to_coords(fine, i, idx)
274 x = idx/2
275 mod2 = mod(idx, 2)
276
277 pt = sum(abs(mod2))
278
279 select case (pt)
280 case (0)
281 i1 = i1 + 1
282 tt%to_fine1(1, i1) = mesh_local_index_from_coords(coarse, x)
283
284 case (1)
285 i2 = i2 + 1
286 tt%to_fine2(1, i2) = mesh_local_index_from_coords(coarse, x)
287 tt%to_fine2(2, i2) = mesh_local_index_from_coords(coarse, x + mod2)
288
289 case (2)
290 i4 = i4 + 1
291 tt%to_fine4(1, i4) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) + mod2(3)])
292 tt%to_fine4(2, i4) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) + mod2(3)])
293 tt%to_fine4(3, i4) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) + mod2(2), x(3) ])
294 tt%to_fine4(4, i4) = mesh_local_index_from_coords(coarse, x + mod2)
295
296 case (3)
297 i8 = i8 + 1
298 tt%to_fine8(1, i8) = mesh_local_index_from_coords(coarse, x)
299 tt%to_fine8(2, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) ])
300 tt%to_fine8(3, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) ])
301 tt%to_fine8(4, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) , x(3) + mod2(3)])
302 tt%to_fine8(5, i8) = mesh_local_index_from_coords(coarse, [x(1) , x(2) + mod2(2), x(3) + mod2(3)])
303 tt%to_fine8(6, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) , x(3) + mod2(3)])
304 tt%to_fine8(7, i8) = mesh_local_index_from_coords(coarse, [x(1) + mod2(1), x(2) + mod2(2), x(3) ])
305 tt%to_fine8(8, i8) = mesh_local_index_from_coords(coarse, x + mod2)
306
307 end select
308
309 end do
310
311 assert(i1 == tt%n_fine1 .and. i2 == tt%n_fine2 .and. i4 == tt%n_fine4 .and. i8 == tt%n_fine8)
312
313
315 end subroutine multigrid_get_transfer_tables
316
317 !---------------------------------------------------------------------------------
320 !---------------------------------------------------------------------------------
321 subroutine multigrid_mesh_half(space, namespace, mesh_in, mesh_out, stencil)
322 class(space_t), intent(in) :: space
323 type(namespace_t), intent(in) :: namespace
324 type(mesh_t), target, intent(in) :: mesh_in
325 type(mesh_t), intent(inout) :: mesh_out
326 type(stencil_t), intent(in) :: stencil
327
328 integer :: idim
329
330 push_sub(multigrid_mesh_half)
331
332 mesh_out%box => mesh_in%box
333 mesh_out%idx%dim = mesh_in%idx%dim
334 mesh_out%use_curvilinear = mesh_in%use_curvilinear
335 mesh_out%masked_periodic_boundaries = mesh_in%masked_periodic_boundaries
336 mesh_out%coord_system => mesh_in%coord_system
337
338 safe_allocate(mesh_out%spacing(1:space%dim))
339 mesh_out%spacing(:) = 2*mesh_in%spacing(:)
340
341 call index_init(mesh_out%idx, space%dim)
342 mesh_out%idx%enlarge(:) = mesh_in%idx%enlarge(:)
343 mesh_out%idx%nr(1,:) = (mesh_in%idx%nr(1,:)+mesh_in%idx%enlarge(:))/2
344 mesh_out%idx%nr(2,:) = (mesh_in%idx%nr(2,:)-mesh_in%idx%enlarge(:))/2
345 mesh_out%idx%ll(:) = mesh_out%idx%nr(2, :) - mesh_out%idx%nr(1, :) + 1
346
347 mesh_out%idx%stride(1) = 1
348 do idim = 2, mesh_out%idx%dim + 1
349 mesh_out%idx%stride(idim) = mesh_out%idx%stride(idim-1) * &
350 (mesh_out%idx%ll(idim-1) + 2*mesh_out%idx%enlarge(idim-1))
351 end do
352
353 call mesh_init_stage_2(mesh_out, namespace, space, mesh_out%box, stencil)
354
355 pop_sub(multigrid_mesh_half)
356 end subroutine multigrid_mesh_half
357
358 !---------------------------------------------------------------------------------
359 subroutine multigrid_mesh_double(space, namespace, mesh_in, mesh_out, stencil)
360 class(space_t), intent(in) :: space
361 type(namespace_t), intent(in) :: namespace
362 type(mesh_t), target, intent(in) :: mesh_in
363 type(mesh_t), intent(inout) :: mesh_out
364 type(stencil_t), intent(in) :: stencil
365
366 integer :: idim
367 push_sub(multigrid_mesh_double)
368
369 mesh_out%box => mesh_in%box
370 mesh_out%idx%dim = mesh_in%idx%dim
371 mesh_out%use_curvilinear = mesh_in%use_curvilinear
372 mesh_out%masked_periodic_boundaries = mesh_in%masked_periodic_boundaries
373 mesh_out%coord_system => mesh_in%coord_system
374
375 safe_allocate(mesh_out%spacing(1:space%dim))
376 mesh_out%spacing(:) = m_half*mesh_in%spacing(:)
377
378 call index_init(mesh_out%idx, space%dim)
379 mesh_out%idx%enlarge(:) = mesh_in%idx%enlarge(:)
380 mesh_out%idx%nr(1,:) = (mesh_in%idx%nr(1,:)+mesh_in%idx%enlarge(:))*2
381 mesh_out%idx%nr(2,:) = (mesh_in%idx%nr(2,:)-mesh_in%idx%enlarge(:))*2
382 ! We need to make the possible number of grid points larger by one for
383 ! the non-periodic dimensions because the spacing is only half of the
384 ! original mesh and thus we could get points in the new boundary
385 ! that are still inside the simulation box.
386 ! For the periodic dimensions, we are anyway commensurate with the size
387 ! of the box, so we are still commensurate when taking twice the number
388 ! of points.
389 do idim = space%periodic_dim + 1, space%dim
390 mesh_out%idx%nr(1, idim) = mesh_out%idx%nr(1, idim) - 1
391 mesh_out%idx%nr(2, idim) = mesh_out%idx%nr(2, idim) + 1
392 end do
393 mesh_out%idx%ll(:) = mesh_out%idx%nr(2, :) - mesh_out%idx%nr(1, :) + 1
394
395 mesh_out%idx%stride(1) = 1
396 do idim = 2, space%dim+1
397 mesh_out%idx%stride(idim) = mesh_out%idx%stride(idim-1) * &
398 (mesh_out%idx%ll(idim-1) + 2*mesh_out%idx%enlarge(idim-1))
399 end do
400
401 call mesh_init_stage_2(mesh_out, namespace, space, mesh_out%box, stencil)
402
403 pop_sub(multigrid_mesh_double)
404 end subroutine multigrid_mesh_double
405
406 ! ---------------------------------------------------------
407 subroutine multigrid_end(mgrid)
408 type(multigrid_t), target, intent(inout) :: mgrid
409
410 integer :: i
411 type(multigrid_level_t), pointer :: level
412
413 push_sub(multigrid_end)
415 safe_deallocate_a(mgrid%sp)
416 safe_deallocate_a(mgrid%ep)
417 safe_deallocate_a(mgrid%ep_part)
418
419 safe_deallocate_a(mgrid%level(0)%tt%fine_i)
420
421 do i = 1, mgrid%n_levels
422 level => mgrid%level(i)
423
424 call derivatives_end(level%der)
425 call mesh_end(level%mesh)
426 safe_deallocate_p(level%mesh)
427 safe_deallocate_p(level%der)
428
429 safe_deallocate_a(level%tt%to_coarse)
430 safe_deallocate_a(level%tt%to_fine1)
431 safe_deallocate_a(level%tt%to_fine2)
432 safe_deallocate_a(level%tt%to_fine4)
433 safe_deallocate_a(level%tt%to_fine8)
434 safe_deallocate_a(level%tt%fine_i)
435 end do
436
437 safe_deallocate_a(mgrid%level)
438
439 pop_sub(multigrid_end)
440 end subroutine multigrid_end
441
442 !---------------------------------------------------------------------------------
443 integer function multigrid_number_of_levels(base_der) result(number)
444 type(derivatives_t), target, intent(in) :: base_der
445
446 type(derivatives_t), pointer :: next_der
447
448 next_der => base_der%coarser
449
450 number = 0
451 do
452 number = number + 1
453 next_der => next_der%coarser
454 if (.not. associated(next_der)) exit
455 end do
456
457 end function multigrid_number_of_levels
458
459
460#include "undef.F90"
461#include "real.F90"
462#include "multigrid_inc.F90"
463
464#include "undef.F90"
465#include "complex.F90"
466#include "multigrid_inc.F90"
467
468end module multigrid_oct_m
469
470!! Local Variables:
471!! mode: f90
472!! coding: utf-8
473!! 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:193
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:534
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_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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:453
subroutine, public zmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, order, set_bc)
Definition: multigrid.F90:1141
subroutine, public dmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
Definition: multigrid.F90:928
subroutine, public zmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb, order)
Definition: multigrid.F90:1343
integer function, public multigrid_number_of_levels(base_der)
Definition: multigrid.F90:537
subroutine, public multigrid_end(mgrid)
Definition: multigrid.F90:501
subroutine, public dmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
Definition: multigrid.F90:703
integer, parameter, public fullweight
Definition: multigrid.F90:157
subroutine, public zmultigrid_fine2coarse_batch(tt, fine_der, coarse_mesh, fineb, coarseb, method_p)
Definition: multigrid.F90:1447
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:415
subroutine, public multigrid_init(mgrid, namespace, space, mesh, der, stencil, mc, nlevels)
Definition: multigrid.F90:184
subroutine, public zmultigrid_fine2coarse(tt, fine_der, coarse_mesh, f_fine, f_coarse, method_p)
Definition: multigrid.F90:1222
subroutine, public dmultigrid_coarse2fine_batch(tt, coarse_der, fine_mesh, coarseb, fineb, order)
Definition: multigrid.F90:824
subroutine, public dmultigrid_coarse2fine(tt, coarse_der, fine_mesh, f_coarse, f_fine, order, set_bc)
Definition: multigrid.F90:622
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:301
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