Afivo
0.3
|
This module contains the geometric multigrid routines that come with Afivo. More...
Data Types | |
interface | mg3_box_op |
Subroutine that performs A * cc(..., i_in) = cc(..., i_out) More... | |
type | mg3_t |
Type to store multigrid options in. More... | |
Functions/Subroutines | |
subroutine, public | mg3_init_mg (mg) |
Check multigrid options or set them to default. More... | |
subroutine, public | mg3_fas_fmg (tree, mg, set_residual, have_guess) |
Perform FAS-FMG cycle (full approximation scheme, full multigrid). Note that this routine needs valid ghost cells (for i_phi) on input, and gives back valid ghost cells on output. More... | |
subroutine, public | mg3_fas_vcycle (tree, mg, set_residual, highest_lvl) |
Perform FAS V-cycle (full approximation scheme). Note that this routine needs valid ghost cells (for i_phi) on input, and gives back valid ghost cells on output. More... | |
subroutine, public | mg3_sides_rb (boxes, id, nb, iv) |
Fill ghost cells near refinement boundaries which preserves diffusive fluxes. More... | |
subroutine | mg3_sides_rb_old (boxes, id, nb, iv) |
Fill ghost cells near refinement boundaries which preserves diffusive fluxes. This routine does not do interpolation of coarse grid values. Basically, we extrapolate from the fine cells to a corner point, and then take the average between this corner point and a coarse neighbor to fill ghost cells for the fine cells. More... | |
subroutine | correct_children (boxes, ids, mg) |
subroutine | update_coarse (tree, lvl, mg) |
subroutine | set_coarse_phi_rhs (tree, lvl, mg) |
This routine performs the same as update_coarse, but it ignores the tmp variable. More... | |
subroutine | init_phi_rhs (tree, mg) |
Set the initial guess for phi and restrict the rhs. More... | |
subroutine | residual_box (box, mg) |
subroutine, public | mg3_auto_gsrb (box, redblack_cntr, mg) |
Based on the box type, apply a Gauss-Seidel relaxation scheme. More... | |
subroutine, public | mg3_auto_op (box, i_out, mg) |
Based on the box type, apply the approriate Laplace operator. More... | |
subroutine, public | mg3_auto_rstr (box_c, box_p, iv, mg) |
Based on the box type, apply the approriate Laplace operator. More... | |
subroutine, public | mg3_auto_corr (box_p, box_c, mg) |
Based on the box type, correct the solution of the children. More... | |
subroutine, public | mg3_set_box_tag (box, mg) |
subroutine, public | mg3_box_corr_lpl (box_p, box_c, mg) |
subroutine, public | mg3_box_gsrb_lpl (box, redblack_cntr, mg) |
Perform Gauss-Seidel relaxation on box for a Laplacian operator. More... | |
subroutine, public | mg3_box_lpl (box, i_out, mg) |
Perform Laplacian operator on a box. More... | |
subroutine, public | mg3_box_rstr_lpl (box_c, box_p, iv, mg) |
Restriction of child box (box_c) to its parent (box_p) More... | |
subroutine, public | mg3_box_gsrb_lpld (box, redblack_cntr, mg) |
Perform Gauss-Seidel relaxation on a box. Epsilon can have a jump at cell faces. More... | |
subroutine, public | mg3_box_lpld (box, i_out, mg) |
Perform Laplacian operator on a box where epsilon varies on cell faces. More... | |
subroutine, public | mg3_box_corr_lpld (box_p, box_c, mg) |
Correct fine grid values based on the change in the coarse grid, in the case of a jump in epsilon. More... | |
subroutine, public | mg3_box_corr_lpllsf (box_p, box_c, mg) |
subroutine, public | mg3_box_gsrb_lpllsf (box, redblack_cntr, mg) |
Perform Gauss-Seidel relaxation on box for a Laplacian operator. More... | |
subroutine, public | mg3_box_lpllsf (box, i_out, mg) |
Perform Laplacian operator on a box. More... | |
subroutine, public | mg3_box_rstr_lpllsf (box_c, box_p, iv, mg) |
Restriction of child box (box_c) to its parent (box_p) More... | |
Variables | |
integer, parameter, public | mg_normal_box = 1 |
Normal box. More... | |
integer, parameter, public | mg_lsf_box = 2 |
Box with an internal boundary. More... | |
integer, parameter, public | mg_ceps_box = 3 |
Box with constant eps /= 1. More... | |
integer, parameter, public | mg_veps_box = 4 |
Box with varying eps (on face) More... | |
This module contains the geometric multigrid routines that come with Afivo.
subroutine, public m_a3_multigrid::mg3_init_mg | ( | type(mg3_t), intent(inout) | mg | ) |
Check multigrid options or set them to default.
[in,out] | mg | Multigrid options |
Definition at line 128 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_fas_fmg | ( | type(a3_t), intent(inout) | tree, |
type(mg3_t), intent(in) | mg, | ||
logical, intent(in) | set_residual, | ||
logical, intent(in) | have_guess | ||
) |
Perform FAS-FMG cycle (full approximation scheme, full multigrid). Note that this routine needs valid ghost cells (for i_phi) on input, and gives back valid ghost cells on output.
[in,out] | tree | Tree to do multigrid on |
[in] | mg | Multigrid options |
[in] | set_residual | If true, store residual in i_tmp |
[in] | have_guess | If false, start from phi = 0 |
Definition at line 162 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_fas_vcycle | ( | type(a3_t), intent(inout) | tree, |
type(mg3_t), intent(in) | mg, | ||
logical, intent(in) | set_residual, | ||
integer, intent(in), optional | highest_lvl | ||
) |
Perform FAS V-cycle (full approximation scheme). Note that this routine needs valid ghost cells (for i_phi) on input, and gives back valid ghost cells on output.
[in,out] | tree | Tree to do multigrid on |
[in] | mg | Multigrid options |
[in] | set_residual | If true, store residual in i_tmp |
[in] | highest_lvl | Maximum level for V-cycle |
Definition at line 208 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_sides_rb | ( | type(box3_t), dimension(:), intent(inout) | boxes, |
integer, intent(in) | id, | ||
integer, intent(in) | nb, | ||
integer, intent(in) | iv | ||
) |
Fill ghost cells near refinement boundaries which preserves diffusive fluxes.
[in,out] | boxes | List of all boxes |
[in] | id | Id of box |
[in] | nb | Ghost cell direction |
[in] | iv | Ghost cell variable |
Definition at line 283 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::mg3_sides_rb_old | ( | type(box3_t), dimension(:), intent(inout) | boxes, |
integer, intent(in) | id, | ||
integer, intent(in) | nb, | ||
integer, intent(in) | iv | ||
) |
Fill ghost cells near refinement boundaries which preserves diffusive fluxes. This routine does not do interpolation of coarse grid values. Basically, we extrapolate from the fine cells to a corner point, and then take the average between this corner point and a coarse neighbor to fill ghost cells for the fine cells.
[in,out] | boxes | List of all boxes |
[in] | id | Id of box |
[in] | nb | Ghost cell direction |
[in] | iv | Ghost cell variable |
Definition at line 391 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::correct_children | ( | type(box3_t), dimension(:), intent(inout) | boxes, |
integer, dimension(:), intent(in) | ids, | ||
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | boxes | List of all boxes |
[in] | ids | Operate on these boxes |
[in] | mg | Multigrid options |
Definition at line 494 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::update_coarse | ( | type(a3_t), intent(inout) | tree, |
integer, intent(in) | lvl, | ||
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | tree | Tree containing full grid |
[in] | lvl | Update coarse values at lvl-1 |
[in] | mg | Multigrid options |
Definition at line 561 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::set_coarse_phi_rhs | ( | type(a3_t), intent(inout) | tree, |
integer, intent(in) | lvl, | ||
type(mg3_t), intent(in) | mg | ||
) |
This routine performs the same as update_coarse, but it ignores the tmp variable.
[in,out] | tree | Tree containing full grid |
[in] | lvl | Update coarse values at lvl-1 |
[in] | mg | Multigrid options |
Definition at line 614 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::init_phi_rhs | ( | type(a3_t), intent(inout) | tree, |
type(mg3_t), intent(in) | mg | ||
) |
Set the initial guess for phi and restrict the rhs.
[in,out] | tree | Full grid |
[in] | mg | Multigrid options |
Definition at line 654 of file m_a3_multigrid.f90.
subroutine m_a3_multigrid::residual_box | ( | type(box3_t), intent(inout) | box, |
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | box | Operate on this box |
[in] | mg | Multigrid options |
Definition at line 679 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_auto_gsrb | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | redblack_cntr, | ||
type(mg3_t), intent(in) | mg | ||
) |
Based on the box type, apply a Gauss-Seidel relaxation scheme.
[in,out] | box | Box to operate on |
[in] | redblack_cntr | Iteration count |
[in] | mg | Multigrid options |
Definition at line 691 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_auto_op | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | i_out, | ||
type(mg3_t), intent(in) | mg | ||
) |
Based on the box type, apply the approriate Laplace operator.
[in,out] | box | Operate on this box |
[in] | i_out | Index of output variable |
[in] | mg | Multigrid options |
Definition at line 709 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_auto_rstr | ( | type(box3_t), intent(in) | box_c, |
type(box3_t), intent(inout) | box_p, | ||
integer, intent(in) | iv, | ||
type(mg3_t), intent(in) | mg | ||
) |
Based on the box type, apply the approriate Laplace operator.
[in] | box_c | Child box |
[in,out] | box_p | Parent box |
[in] | iv | Index of variable |
[in] | mg | Multigrid options |
Definition at line 727 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_auto_corr | ( | type(box3_t), intent(in) | box_p, |
type(box3_t), intent(inout) | box_c, | ||
type(mg3_t), intent(in) | mg | ||
) |
Based on the box type, correct the solution of the children.
[in,out] | box_c | Child box |
[in] | box_p | Parent box |
[in] | mg | Multigrid options |
Definition at line 745 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_set_box_tag | ( | type(box3_t), intent(inout) | box, |
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | box | Box to operate on |
[in] | mg | Multigrid options |
Definition at line 762 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_corr_lpl | ( | type(box3_t), intent(in) | box_p, |
type(box3_t), intent(inout) | box_c, | ||
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | box_c | Child box |
[in] | box_p | Parent box |
[in] | mg | Multigrid options |
Definition at line 793 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_gsrb_lpl | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | redblack_cntr, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Gauss-Seidel relaxation on box for a Laplacian operator.
[in,out] | box | Box to operate on |
[in] | redblack_cntr | Iteration counter |
[in] | mg | Multigrid options |
Definition at line 803 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_lpl | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | i_out, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Laplacian operator on a box.
[in,out] | box | Box to operate on |
[in] | i_out | Index of variable to store Laplacian in |
[in] | mg | Multigrid options |
Definition at line 834 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_rstr_lpl | ( | type(box3_t), intent(in) | box_c, |
type(box3_t), intent(inout) | box_p, | ||
integer, intent(in) | iv, | ||
type(mg3_t), intent(in) | mg | ||
) |
Restriction of child box (box_c) to its parent (box_p)
[in] | box_c | Child box to restrict |
[in,out] | box_p | Parent box to restrict to |
[in] | iv | Variable to restrict |
[in] | mg | Multigrid options |
Definition at line 859 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_gsrb_lpld | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | redblack_cntr, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Gauss-Seidel relaxation on a box. Epsilon can have a jump at cell faces.
[in,out] | box | Box to operate on |
[in] | redblack_cntr | Iteration counter |
[in] | mg | Multigrid options |
Definition at line 873 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_lpld | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | i_out, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Laplacian operator on a box where epsilon varies on cell faces.
[in,out] | box | Box to operate on |
[in] | i_out | Index of variable to store Laplacian in |
[in] | mg | Multigrid options |
Definition at line 910 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_corr_lpld | ( | type(box3_t), intent(in) | box_p, |
type(box3_t), intent(inout) | box_c, | ||
type(mg3_t), intent(in) | mg | ||
) |
Correct fine grid values based on the change in the coarse grid, in the case of a jump in epsilon.
[in,out] | box_c | Child box |
[in] | box_p | Parent box |
[in] | mg | Multigrid options |
Definition at line 945 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_corr_lpllsf | ( | type(box3_t), intent(in) | box_p, |
type(box3_t), intent(inout) | box_c, | ||
type(mg3_t), intent(in) | mg | ||
) |
[in,out] | box_c | Child box |
[in] | box_p | Parent box |
[in] | mg | Multigrid options |
Definition at line 1024 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_gsrb_lpllsf | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | redblack_cntr, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Gauss-Seidel relaxation on box for a Laplacian operator.
[in,out] | box | Box to operate on |
[in] | redblack_cntr | Iteration counter |
[in] | mg | Multigrid options |
Definition at line 1077 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_lpllsf | ( | type(box3_t), intent(inout) | box, |
integer, intent(in) | i_out, | ||
type(mg3_t), intent(in) | mg | ||
) |
Perform Laplacian operator on a box.
[in,out] | box | Box to operate on |
[in] | i_out | Index of variable to store Laplacian in |
[in] | mg | Multigrid options |
Definition at line 1129 of file m_a3_multigrid.f90.
subroutine, public m_a3_multigrid::mg3_box_rstr_lpllsf | ( | type(box3_t), intent(in) | box_c, |
type(box3_t), intent(inout) | box_p, | ||
integer, intent(in) | iv, | ||
type(mg3_t), intent(in) | mg | ||
) |
Restriction of child box (box_c) to its parent (box_p)
[in] | box_c | Child box to restrict |
[in,out] | box_p | Parent box to restrict to |
[in] | iv | Variable to restrict |
[in] | mg | Multigrid options |
Definition at line 1176 of file m_a3_multigrid.f90.
integer, parameter, public m_a3_multigrid::mg_normal_box = 1 |
Normal box.
Definition at line 11 of file m_a3_multigrid.f90.
integer, parameter, public m_a3_multigrid::mg_lsf_box = 2 |
Box with an internal boundary.
Definition at line 12 of file m_a3_multigrid.f90.
integer, parameter, public m_a3_multigrid::mg_ceps_box = 3 |
Box with constant eps /= 1.
Definition at line 13 of file m_a3_multigrid.f90.
integer, parameter, public m_a3_multigrid::mg_veps_box = 4 |
Box with varying eps (on face)
Definition at line 14 of file m_a3_multigrid.f90.