Afivo
0.3
|
Module to solve elliptic PDEs on the coarse grid. This module contains an interface to Hypre, assuming Hypre is compiled with OpenMP and without MPI. More...
Data Types | |
interface | HYPRE_StructMatrixSetBoxValues |
Functions/Subroutines | |
subroutine, public | coarse_solver_initialize (tree, mg) |
Initialize the coarse grid solver. More... | |
subroutine | hypre_set_matrix (tree, mg) |
Set matrix type and store coefficients. More... | |
subroutine, public | coarse_solver_update_matrix (tree, mg) |
Update matrix coefficients. More... | |
subroutine, public | coarse_solver_destroy (cs) |
De-allocate storage for all coarse solver components. More... | |
subroutine | hypre_destroy_solver (cs) |
De-allocate storage for solver. More... | |
subroutine | hypre_create_grid (grid, nx, periodic) |
subroutine | hypre_create_vector (grid, nx, vec) |
Create a Hypre vector. More... | |
subroutine, public | coarse_solver_set_rhs_phi (tree, mg) |
Set the right-hand side and copy phi from the tree. Also move the boundary conditions for phi to the rhs. More... | |
subroutine, public | coarse_solver_get_phi (tree, mg) |
Copy solution from coarse solver to the tree. More... | |
subroutine | hypre_create_matrix (A, grid, stencil_size, offsets, symmetric) |
Create a symmetric matrix on a grid. More... | |
subroutine | hypre_prepare_solve (cs) |
subroutine, public | coarse_solver (cs, num_iterations) |
subroutine | stencil_handle_boundaries (box, mg, stencil, bc_to_rhs) |
Incorporate boundary conditions into stencil. More... | |
real(dp) function, dimension(dtimes(box%n_cell)), public | mg_lsf_boundary_value (box, mg) |
Compute boundary value for internal boundaries. More... | |
Variables | |
integer, parameter | hypre_parcsr = 5555 |
integer, parameter, public | coarse_solver_hypre_smg = 1 |
integer, parameter, public | coarse_solver_hypre_pfmg = 2 |
integer, parameter, public | coarse_solver_hypre_cycred = 3 |
integer, parameter | max_stencil_size = 2*NDIM + 1 |
integer, dimension(3, max_stencil_size), parameter | stencil_offsets = reshape([0, -1, 1], [1, max_stencil_size]) |
integer, parameter | mpi_comm_world = 0 |
integer, parameter | num_procs = 1 |
integer, parameter | myid = 0 |
Module to solve elliptic PDEs on the coarse grid. This module contains an interface to Hypre, assuming Hypre is compiled with OpenMP and without MPI.
subroutine, public m_coarse_solver::coarse_solver_initialize | ( | type(af_t), intent(inout) | tree, |
type(mg_t), intent(inout) | mg | ||
) |
Initialize the coarse grid solver.
[in,out] | tree | Tree to do multigrid on |
Definition at line 71 of file m_coarse_solver.f90.
|
private |
Set matrix type and store coefficients.
Definition at line 104 of file m_coarse_solver.f90.
subroutine, public m_coarse_solver::coarse_solver_update_matrix | ( | type(af_t), intent(inout) | tree, |
type(mg_t), intent(inout) | mg | ||
) |
Update matrix coefficients.
Definition at line 197 of file m_coarse_solver.f90.
subroutine, public m_coarse_solver::coarse_solver_destroy | ( | type(coarse_solve_t), intent(inout) | cs | ) |
De-allocate storage for all coarse solver components.
Definition at line 208 of file m_coarse_solver.f90.
|
private |
De-allocate storage for solver.
Definition at line 222 of file m_coarse_solver.f90.
|
private |
[in] | nx | Size of grid |
[in] | periodic | Whether the dimension is periodic |
Definition at line 238 of file m_coarse_solver.f90.
|
private |
Create a Hypre vector.
Definition at line 268 of file m_coarse_solver.f90.
subroutine, public m_coarse_solver::coarse_solver_set_rhs_phi | ( | type(af_t), intent(inout) | tree, |
type(mg_t), intent(in) | mg | ||
) |
Set the right-hand side and copy phi from the tree. Also move the boundary conditions for phi to the rhs.
Definition at line 286 of file m_coarse_solver.f90.
subroutine, public m_coarse_solver::coarse_solver_get_phi | ( | type(af_t), intent(inout) | tree, |
type(mg_t), intent(in) | mg | ||
) |
Copy solution from coarse solver to the tree.
Definition at line 341 of file m_coarse_solver.f90.
|
private |
Create a symmetric matrix on a grid.
Definition at line 361 of file m_coarse_solver.f90.
|
private |
Definition at line 393 of file m_coarse_solver.f90.
subroutine, public m_coarse_solver::coarse_solver | ( | type(coarse_solve_t), intent(in) | cs, |
integer, intent(out) | num_iterations | ||
) |
Definition at line 421 of file m_coarse_solver.f90.
|
private |
Incorporate boundary conditions into stencil.
Definition at line 442 of file m_coarse_solver.f90.
real(dp) function, dimension(dtimes(box%n_cell)), public m_coarse_solver::mg_lsf_boundary_value | ( | type(box_t), intent(in) | box, |
type(mg_t), intent(in) | mg | ||
) |
Compute boundary value for internal boundaries.
Definition at line 494 of file m_coarse_solver.f90.
|
private |
Definition at line 12 of file m_coarse_solver.f90.
integer, parameter, public m_coarse_solver::coarse_solver_hypre_smg = 1 |
Definition at line 15 of file m_coarse_solver.f90.
integer, parameter, public m_coarse_solver::coarse_solver_hypre_pfmg = 2 |
Definition at line 18 of file m_coarse_solver.f90.
integer, parameter, public m_coarse_solver::coarse_solver_hypre_cycred = 3 |
Definition at line 21 of file m_coarse_solver.f90.
|
private |
Definition at line 24 of file m_coarse_solver.f90.
|
private |
Definition at line 28 of file m_coarse_solver.f90.
|
private |
Definition at line 42 of file m_coarse_solver.f90.
|
private |
Definition at line 43 of file m_coarse_solver.f90.
|
private |
Definition at line 44 of file m_coarse_solver.f90.