Afivo  0.3
m_coarse_solver.f90
1 #include "../src/cpp_macros.h"
2 
5  use m_af_types
6  use m_af_stencil
7 
8  implicit none
9  private
10 
11  ! Hypre CSR format
12  integer, parameter :: hypre_parcsr = 5555
13 
14  ! Semi-coarsening multigrid (more robust and expensive)
15  integer, parameter, public :: coarse_solver_hypre_smg = 1
16 
17  ! Multigrid with point-wise smoother (fast/cheap)
18  integer, parameter, public :: coarse_solver_hypre_pfmg = 2
19 
20  ! Cyclic reduction solver
21  integer, parameter, public :: coarse_solver_hypre_cycred = 3
22 
23  ! Size of the stencil (only direct neighbors)
24  integer, parameter :: max_stencil_size = 2*ndim + 1
25 
26  ! Offsets for the stencil elements
27 #if NDIM == 1
28  integer, parameter :: stencil_offsets(1, max_stencil_size) = reshape([0, &
29  -1, 1], [1, max_stencil_size])
30 #elif NDIM == 2
31  integer, parameter :: stencil_offsets(2, max_stencil_size) = reshape([0, 0, &
32  -1, 0, 1, 0, &
33  0, -1, 0, 1], [2, max_stencil_size])
34 #elif NDIM == 3
35  integer, parameter :: stencil_offsets(3, max_stencil_size) = reshape([0, 0, 0, &
36  -1, 0, 0, 1, 0, 0, &
37  0, -1, 0, 0, 1, 0, &
38  0, 0, -1, 0, 0, 1], [3, max_stencil_size])
39 #endif
40 
41  ! MPI stubs
42  integer, parameter :: mpi_comm_world = 0
43  integer, parameter :: num_procs = 1
44  integer, parameter :: myid = 0
45 
46  interface
47  subroutine hypre_structmatrixsetboxvalues(matrix, ilower, iupper, nentries, &
48  entries, values, ierr)
49  import
50  type(c_ptr), intent(in) :: matrix
51  integer, intent(in) :: ilower(*), iupper(*)
52  integer, intent(in) :: nentries, entries(*)
53  real(dp), intent(in) :: values(*)
54  integer, intent(out) :: ierr
55  end subroutine hypre_structmatrixsetboxvalues
56  end interface
57 
58  public :: coarse_solver_initialize
59  public :: coarse_solver_destroy
60  public :: coarse_solver_set_rhs_phi
61  public :: coarse_solver_get_phi
62  public :: coarse_solver
63  public :: coarse_solver_update_matrix
64 
65  ! Could be moved to another module at some point
66  public :: mg_lsf_boundary_value
67 
68 contains
69 
71  subroutine coarse_solver_initialize(tree, mg)
72  type(af_t), intent(inout) :: tree
73  type(mg_t), intent(inout) :: mg
74  integer :: nx(ndim), ierr
75 
76  nx = tree%coarse_grid_size(1:ndim)
77 
78  if (.not. tree%ready) error stop "coarse_solver_initialize: tree not ready"
79  if (any(nx == -1)) &
80  error stop "coarse_solver_initialize: coarse_grid_size not set"
81 
82  call hypre_initialize(ierr)
83 
84  ! Construct grid and vectors
85  call hypre_create_grid(mg%csolver%grid, nx, tree%periodic)
86  call hypre_create_vector(mg%csolver%grid, nx, mg%csolver%phi)
87  call hypre_create_vector(mg%csolver%grid, nx, mg%csolver%rhs)
88  call hypre_set_matrix(tree, mg)
89 
90  if (mg%csolver%solver_type <= 0) then
91  if (ndim == 1) then
92  ! PFMG cannot be used in 1D
93  mg%csolver%solver_type = coarse_solver_hypre_cycred
94  else
95  mg%csolver%solver_type = coarse_solver_hypre_pfmg
96  end if
97  end if
98 
99  call hypre_prepare_solve(mg%csolver)
100 
101  end subroutine coarse_solver_initialize
102 
104  subroutine hypre_set_matrix(tree, mg)
105  type(af_t), intent(inout) :: tree
106  type(mg_t), intent(inout) :: mg
107  integer :: n, n_boxes, nc, id, IJK
108  integer :: cnt, stencil_size, ierr
109  integer :: ilo(NDIM), ihi(NDIM)
110  real(dp), allocatable :: coeffs(:, :)
111  real(dp) :: full_coeffs(2*NDIM+1, DTIMES(tree%n_cell))
112  integer, allocatable :: stencil_ix(:)
113  integer, parameter :: zero_to_n(max_stencil_size) = [(i, i=0, max_stencil_size-1)]
114 
115  nc = tree%n_cell
116  n_boxes = size(tree%lvls(1)%ids)
117 
118  if (tree%coord_t == af_cyl .or. mg%i_lsf /= -1) then
119  ! The symmetry option does not seem to work well with axisymmetric
120  ! problems. It also doesn't work with a level set function internal
121  ! boundary.
122  mg%csolver%symmetric = 0
123  end if
124 
125  if (mg%csolver%symmetric == 1) then
126  stencil_size = ndim + 1
127  allocate(stencil_ix(stencil_size))
128 #if NDIM == 1
129  stencil_ix = [1, 3]
130 #elif NDIM == 2
131  stencil_ix = [1, 3, 5]
132 #elif NDIM == 3
133  stencil_ix = [1, 3, 5, 7]
134 #endif
135  else
136  stencil_size = 2*ndim + 1
137  allocate(stencil_ix(stencil_size))
138 #if NDIM == 1
139  stencil_ix = [1, 2, 3]
140 #elif NDIM == 2
141  stencil_ix = [1, 2, 3, 4, 5]
142 #elif NDIM == 3
143  stencil_ix = [1, 2, 3, 4, 5, 6, 7]
144 #endif
145  end if
146 
147  ! Construct the matrix
148  call hypre_create_matrix(mg%csolver%matrix, mg%csolver%grid, &
149  stencil_size, stencil_offsets(:, stencil_ix), mg%csolver%symmetric)
150 
151  if (.not. allocated(mg%csolver%bc_to_rhs)) then
152  allocate(mg%csolver%bc_to_rhs(nc**(ndim-1), af_num_neighbors, n_boxes))
153  end if
154  if (.not. allocated(mg%csolver%lsf_fac)) then
155  allocate(mg%csolver%lsf_fac(dtimes(nc), n_boxes))
156  end if
157  allocate(coeffs(stencil_size, tree%n_cell**ndim))
158 
159  mg%csolver%bc_to_rhs = 0.0_dp
160  mg%csolver%lsf_fac = 0.0_dp
161  coeffs = 0.0_dp
162 
163  do n = 1, size(tree%lvls(1)%ids)
164  id = tree%lvls(1)%ids(n)
165 
166  associate(box => tree%boxes(id))
167  call af_stencil_get_box(box, mg%operator_key, full_coeffs)
168  call stencil_handle_boundaries(box, mg, full_coeffs, &
169  mg%csolver%bc_to_rhs(:, :, n))
170 
171  ! This assumes that correction factors due to a level set function are
172  ! stored in the stencil%f array
173  i = af_stencil_index(box, mg%operator_key)
174  if (allocated(box%stencils(i)%f)) then
175  mg%csolver%lsf_fac(dtimes(:), n) = box%stencils(i)%f
176  end if
177  end associate
178 
179  cnt = 0
180  do kji_do(1, nc)
181  cnt = cnt + 1
182  coeffs(:, cnt) = full_coeffs(stencil_ix, ijk)
183  end do; close_do
184 
185  ilo = (tree%boxes(id)%ix - 1) * nc + 1
186  ihi = ilo + nc - 1
187  call hypre_structmatrixsetboxvalues(mg%csolver%matrix, ilo, ihi, &
188  stencil_size, zero_to_n(1:stencil_size), coeffs, ierr)
189  if (ierr /= 0) error stop "HYPRE_StructMatrixSetBoxValues failed"
190  end do
191 
192  call hypre_structmatrixassemble(mg%csolver%matrix, ierr)
193 
194  end subroutine hypre_set_matrix
195 
197  subroutine coarse_solver_update_matrix(tree, mg)
198  type(af_t), intent(inout) :: tree
199  type(mg_t), intent(inout) :: mg
200  integer :: ierr
201 
202  call hypre_structmatrixdestroy(mg%csolver%matrix, ierr)
203  call hypre_set_matrix(tree, mg)
204  call hypre_prepare_solve(mg%csolver)
205  end subroutine coarse_solver_update_matrix
206 
208  subroutine coarse_solver_destroy(cs)
209  type(coarse_solve_t), intent(inout) :: cs
210  integer :: ierr
211 
212  call hypre_structgriddestroy(cs%grid, ierr)
213  call hypre_structmatrixdestroy(cs%matrix, ierr)
214  call hypre_structvectordestroy(cs%rhs, ierr)
215  call hypre_structvectordestroy(cs%phi, ierr)
216 
217  call hypre_destroy_solver(cs)
218  call hypre_finalize(ierr)
219  end subroutine coarse_solver_destroy
220 
222  subroutine hypre_destroy_solver(cs)
223  type(coarse_solve_t), intent(inout) :: cs
224  integer :: ierr
225 
226  select case (cs%solver_type)
227  case (coarse_solver_hypre_cycred)
228  call hypre_structcycreddestroy(cs%solver, ierr)
229  case (coarse_solver_hypre_smg)
230  call hypre_structsmgdestroy(cs%solver, ierr)
231  case (coarse_solver_hypre_pfmg)
232  call hypre_structpfmgdestroy(cs%solver, ierr)
233  case default
234  error stop "hypre_solver_destroy: unknown solver type"
235  end select
236  end subroutine hypre_destroy_solver
237 
238  subroutine hypre_create_grid(grid, nx, periodic)
239  type(c_ptr), intent(out) :: grid
240  integer, intent(in) :: nx(NDIM)
241  logical, intent(in) :: periodic(NDIM)
242  integer :: ierr, ilo(NDIM)
243  integer :: period(NDIM)
244 
245  ilo(:) = 1
246 
247  ! Create an empty grid object
248  call hypre_structgridcreate(mpi_comm_world, ndim, grid, ierr)
249 
250  ! Add a new box to the grid
251  call hypre_structgridsetextents(grid, ilo, nx, ierr)
252 
253  ! Set periodic
254  where (periodic)
255  period = nx
256  elsewhere
257  period = 0
258  end where
259 
260  call hypre_structgridsetperiodic(grid, period, ierr)
261 
262  ! This is a collective call finalizing the grid assembly. The grid is now
263  ! ``ready to be used''
264  call hypre_structgridassemble(grid, ierr)
265  end subroutine hypre_create_grid
266 
268  subroutine hypre_create_vector(grid, nx, vec)
269  type(c_ptr), intent(in) :: grid
270  integer, intent(in) :: nx(:)
271  type(c_ptr), intent(out) :: vec
272  integer :: ierr
273 
274  ! Create an empty vector object */
275  call hypre_structvectorcreate(mpi_comm_world, grid, vec, ierr)
276 
277  ! Indicate that the vector coefficients are ready to be set */
278  call hypre_structvectorinitialize(vec, ierr)
279 
280  ! Finalize the construction of the vector before using
281  call hypre_structvectorassemble(vec, ierr)
282  end subroutine hypre_create_vector
283 
286  subroutine coarse_solver_set_rhs_phi(tree, mg)
287  type(af_t), intent(inout) :: tree
288  type(mg_t), intent(in) :: mg
289  integer :: n, nb, nc, id, bc_type, ierr
290  integer :: ilo(ndim), ihi(ndim), ix
291  real(dp) :: tmp(dtimes(tree%n_cell))
292  real(dp) :: bc_val(tree%n_cell**(ndim-1))
293 
294  nc = tree%n_cell
295 
296  do n = 1, size(tree%lvls(1)%ids)
297  id = tree%lvls(1)%ids(n)
298 
299  tmp = tree%boxes(id)%cc(dtimes(1:nc), mg%i_rhs)
300 
301  ! Add contribution of boundary conditions to rhs
302  do nb = 1, af_num_neighbors
303  if (tree%boxes(id)%neighbors(nb) < af_no_box) then
304  ix = tree%boxes(id)%nb_to_bc_index(nb)
305  call mg%sides_bc(tree%boxes(id), nb, mg%i_phi, &
306  tree%boxes(id)%bc_coords(:, :, ix), &
307  bc_val, bc_type)
308  tree%boxes(id)%bc_val(:, mg%i_phi, ix) = bc_val
309  tree%boxes(id)%bc_type(mg%i_phi, ix) = bc_type
310 
311  ! Get index range near neighbor
312  call af_get_index_bc_inside(nb, nc, 1, ilo, ihi)
313 
314  ! Use the stored arrays mg%csolver%bc_to_rhs to convert the value
315  ! at the boundary to the rhs
316  tmp(dslice(ilo, ihi)) = &
317  tmp(dslice(ilo, ihi)) + &
318  reshape(mg%csolver%bc_to_rhs(:, nb, n) * pack(bc_val, .true.), &
319  [ihi - ilo + 1])
320  end if
321  end do
322 
323  ! Add contribution of level-set function to rhs
324  if (mg%i_lsf /= -1) then
325  tmp = tmp + mg%csolver%lsf_fac(dtimes(:), n) * &
326  mg_lsf_boundary_value(tree%boxes(id), mg)
327  end if
328 
329  ilo = (tree%boxes(id)%ix - 1) * nc + 1
330  ihi = ilo + nc - 1
331  call hypre_structvectorsetboxvalues(mg%csolver%rhs, ilo, ihi, &
332  pack(tmp, .true.), ierr)
333 
334  tmp = tree%boxes(id)%cc(dtimes(1:nc), mg%i_phi)
335  call hypre_structvectorsetboxvalues(mg%csolver%phi, ilo, ihi, &
336  pack(tmp, .true.), ierr)
337  end do
338  end subroutine coarse_solver_set_rhs_phi
339 
341  subroutine coarse_solver_get_phi(tree, mg)
342  type(af_t), intent(inout) :: tree
343  type(mg_t), intent(in) :: mg
344  integer :: n, nc, id, ierr
345  integer :: ilo(ndim), ihi(ndim)
346  real(dp) :: tmp(tree%n_cell**ndim)
347 
348  nc = tree%n_cell
349 
350  do n = 1, size(tree%lvls(1)%ids)
351  id = tree%lvls(1)%ids(n)
352  ilo = (tree%boxes(id)%ix - 1) * nc + 1
353  ihi = ilo + nc - 1
354 
355  call hypre_structvectorgetboxvalues(mg%csolver%phi, ilo, ihi, tmp, ierr)
356  tree%boxes(id)%cc(dtimes(1:nc), mg%i_phi) = reshape(tmp, [dtimes(nc)])
357  end do
358  end subroutine coarse_solver_get_phi
359 
361  subroutine hypre_create_matrix(A, grid, stencil_size, offsets, symmetric)
362  type(c_ptr), intent(out) :: A
363  type(c_ptr), intent(in) :: grid
364  integer, intent(in) :: stencil_size
365  integer, intent(in) :: offsets(NDIM, stencil_size)
366  integer, intent(in) :: symmetric
367  type(c_ptr) :: stencil
368  integer :: i, ierr
369 
370  ! Create an empty stencil object (symmetric)
371  call hypre_structstencilcreate(ndim, stencil_size, stencil, ierr)
372 
373  ! Assign stencil entries */
374  do i = 1, stencil_size
375  call hypre_structstencilsetelement(stencil, i-1, offsets(:, i), ierr)
376  end do
377 
378  ! Create an empty matrix object */
379  call hypre_structmatrixcreate(mpi_comm_world, grid, stencil, a, ierr)
380 
381  ! Use symmetric storage? */
382  call hypre_structmatrixsetsymmetric(a, symmetric, ierr)
383 
384  ! Indicate that the matrix coefficients are ready to be set */
385  call hypre_structmatrixinitialize(a, ierr)
386 
387  call hypre_structstencildestroy(stencil, ierr)
388 
389  end subroutine hypre_create_matrix
390 
391  ! Prepare to solve the system. The coefficient data in b and x is ignored
392  ! here, but information about the layout of the data may be used.
393  subroutine hypre_prepare_solve(cs)
394  type(coarse_solve_t), intent(inout) :: cs
395  integer :: ierr
396 
397  select case (cs%solver_type)
398  case (coarse_solver_hypre_cycred)
399  call hypre_structcycredcreate(mpi_comm_world, cs%solver, ierr)
400  call hypre_structcycredsetup(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
401  case (coarse_solver_hypre_smg)
402  call hypre_structsmgcreate(mpi_comm_world, cs%solver, ierr)
403  call hypre_structsmgsetmaxiter(cs%solver, cs%max_iterations, ierr)
404  call hypre_structsmgsettol(cs%solver, cs%tolerance, ierr)
405  call hypre_structsmgsetnumprerelax(cs%solver, cs%n_cycle_down, ierr)
406  call hypre_structsmgsetnumpostrelax(cs%solver, cs%n_cycle_up, ierr)
407  call hypre_structsmgsetup(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
408  case (coarse_solver_hypre_pfmg)
409  call hypre_structpfmgcreate(mpi_comm_world, cs%solver, ierr)
410  call hypre_structpfmgsetmaxiter(cs%solver, cs%max_iterations, ierr)
411  call hypre_structpfmgsettol(cs%solver, cs%tolerance, ierr)
412  call hypre_structpfmgsetnumprerelax(cs%solver, cs%n_cycle_down, ierr)
413  call hypre_structpfmgsetnumpostrelax(cs%solver, cs%n_cycle_up, ierr)
414  call hypre_structpfmgsetup(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
415  case default
416  error stop "hypre_prepare_solve: unknown solver type"
417  end select
418  end subroutine hypre_prepare_solve
419 
420  ! Solve the system A x = b
421  subroutine coarse_solver(cs, num_iterations)
422  type(coarse_solve_t), intent(in) :: cs
423  integer, intent(out) :: num_iterations
424  integer :: ierr
425 
426  select case (cs%solver_type)
427  case (coarse_solver_hypre_cycred)
428  call hypre_structcycredsolve(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
429  ! call HYPRE_StructCycRedGetNumIterations(cs%solver, num_iterations, ierr)
430  case (coarse_solver_hypre_smg)
431  call hypre_structsmgsolve(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
432  call hypre_structsmggetnumiterations(cs%solver, num_iterations, ierr)
433  case (coarse_solver_hypre_pfmg)
434  call hypre_structpfmgsolve(cs%solver, cs%matrix, cs%rhs, cs%phi, ierr)
435  call hypre_structpfmggetnumiteration(cs%solver, num_iterations, ierr)
436  case default
437  error stop "coarse_solver: unknown solver type"
438  end select
439  end subroutine coarse_solver
440 
442  subroutine stencil_handle_boundaries(box, mg, stencil, bc_to_rhs)
443  type(box_t), intent(in) :: box
444  type(mg_t), intent(in) :: mg
445  real(dp), intent(inout) :: stencil(2*NDIM+1, DTIMES(box%n_cell))
446  real(dp), intent(inout) :: bc_to_rhs(box%n_cell**(NDIM-1), af_num_neighbors)
447  integer :: nb, nc, lo(NDIM), hi(NDIM)
448  integer :: nb_id, nb_dim, bc_type
449  real(dp) :: coords(NDIM, box%n_cell**(NDIM-1))
450  real(dp) :: bc_val(box%n_cell**(NDIM-1))
451 
452  bc_to_rhs = 0.0_dp
453  nc = box%n_cell
454 
455  do nb = 1, af_num_neighbors
456  nb_id = box%neighbors(nb)
457 
458  if (nb_id < af_no_box) then
459  nb_dim = af_neighb_dim(nb)
460  call af_get_face_coords(box, nb, coords)
461  call mg%sides_bc(box, nb, mg%i_phi, coords, bc_val, bc_type)
462 
463  ! Determine index range next to boundary
464  call af_get_index_bc_inside(nb, nc, 1, lo, hi)
465 
466  select case (bc_type)
467  case (af_bc_dirichlet)
468  ! Dirichlet value at cell face, so compute gradient over h/2
469  ! E.g. 1 -2 1 becomes 0 -3 1 for a 1D Laplacian
470  ! The boundary condition is incorporated in the right-hand side
471  stencil(1, dslice(lo, hi)) = &
472  stencil(1, dslice(lo, hi)) - &
473  stencil(nb+1, dslice(lo, hi))
474  bc_to_rhs(:, nb) = pack(-2 * stencil(nb+1, dslice(lo, hi)), .true.)
475  stencil(nb+1, dslice(lo, hi)) = 0.0_dp
476  case (af_bc_neumann)
477  ! E.g. 1 -2 1 becomes 0 -1 1 for a 1D Laplacian
478  stencil(1, dslice(lo, hi)) = &
479  stencil(1, dslice(lo, hi)) + &
480  stencil(nb+1, dslice(lo, hi))
481  bc_to_rhs(:, nb) = &
482  -pack(stencil(nb+1, dslice(lo, hi)) * &
483  box%dr(nb_dim), .true.) * af_neighb_high_pm(nb)
484  stencil(nb+1, dslice(lo, hi)) = 0.0_dp
485  case default
486  error stop "mg_box_lpl_stencil: unsupported boundary condition"
487  end select
488  end if
489  end do
490 
491  end subroutine stencil_handle_boundaries
492 
494  function mg_lsf_boundary_value(box, mg) result(bc)
495  type(box_t), intent(in) :: box
496  type(mg_t), intent(in) :: mg
497  real(dp) :: bc(dtimes(box%n_cell))
498  integer :: ijk, nc
499  real(dp) :: r(ndim)
500 
501  nc = box%n_cell
502 
503  if (associated(mg%lsf_boundary_function)) then
504  do kji_do(1,nc)
505  r = af_r_cc(box, [ijk])
506  bc(ijk) = mg%lsf_boundary_function(r)
507  end do; close_do
508  else
509  bc = mg%lsf_boundary_value
510  end if
511  end function mg_lsf_boundary_value
512 
513 end module m_coarse_solver
This module contains functionality for dealing with numerical stencils.
Definition: m_af_stencil.f90:3
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Definition: m_af_types.f90:4
Module to solve elliptic PDEs on the coarse grid. This module contains an interface to Hypre,...
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
Definition: m_af_types.f90:326
Type to store multigrid options in.
Definition: m_af_types.f90:572