Afivo  0.3
Data Types | Functions/Subroutines | Variables
m_a3_multigrid Module Reference

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...
 

Detailed Description

This module contains the geometric multigrid routines that come with Afivo.

Function/Subroutine Documentation

§ mg3_init_mg()

subroutine, public m_a3_multigrid::mg3_init_mg ( type(mg3_t), intent(inout)  mg)

Check multigrid options or set them to default.

Parameters
[in,out]mgMultigrid options

Definition at line 128 of file m_a3_multigrid.f90.

§ mg3_fas_fmg()

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.

Parameters
[in,out]treeTree to do multigrid on
[in]mgMultigrid options
[in]set_residualIf true, store residual in i_tmp
[in]have_guessIf false, start from phi = 0

Definition at line 162 of file m_a3_multigrid.f90.

§ mg3_fas_vcycle()

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.

Parameters
[in,out]treeTree to do multigrid on
[in]mgMultigrid options
[in]set_residualIf true, store residual in i_tmp
[in]highest_lvlMaximum level for V-cycle

Definition at line 208 of file m_a3_multigrid.f90.

§ mg3_sides_rb()

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.

Parameters
[in,out]boxesList of all boxes
[in]idId of box
[in]nbGhost cell direction
[in]ivGhost cell variable

Definition at line 283 of file m_a3_multigrid.f90.

§ mg3_sides_rb_old()

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.

Parameters
[in,out]boxesList of all boxes
[in]idId of box
[in]nbGhost cell direction
[in]ivGhost cell variable

Definition at line 391 of file m_a3_multigrid.f90.

§ correct_children()

subroutine m_a3_multigrid::correct_children ( type(box3_t), dimension(:), intent(inout)  boxes,
integer, dimension(:), intent(in)  ids,
type(mg3_t), intent(in)  mg 
)
Parameters
[in,out]boxesList of all boxes
[in]idsOperate on these boxes
[in]mgMultigrid options

Definition at line 494 of file m_a3_multigrid.f90.

§ update_coarse()

subroutine m_a3_multigrid::update_coarse ( type(a3_t), intent(inout)  tree,
integer, intent(in)  lvl,
type(mg3_t), intent(in)  mg 
)
Parameters
[in,out]treeTree containing full grid
[in]lvlUpdate coarse values at lvl-1
[in]mgMultigrid options

Definition at line 561 of file m_a3_multigrid.f90.

§ set_coarse_phi_rhs()

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.

Parameters
[in,out]treeTree containing full grid
[in]lvlUpdate coarse values at lvl-1
[in]mgMultigrid options

Definition at line 614 of file m_a3_multigrid.f90.

§ init_phi_rhs()

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.

Parameters
[in,out]treeFull grid
[in]mgMultigrid options

Definition at line 654 of file m_a3_multigrid.f90.

§ residual_box()

subroutine m_a3_multigrid::residual_box ( type(box3_t), intent(inout)  box,
type(mg3_t), intent(in)  mg 
)
Parameters
[in,out]boxOperate on this box
[in]mgMultigrid options

Definition at line 679 of file m_a3_multigrid.f90.

§ mg3_auto_gsrb()

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.

Parameters
[in,out]boxBox to operate on
[in]redblack_cntrIteration count
[in]mgMultigrid options

Definition at line 691 of file m_a3_multigrid.f90.

§ mg3_auto_op()

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.

Parameters
[in,out]boxOperate on this box
[in]i_outIndex of output variable
[in]mgMultigrid options

Definition at line 709 of file m_a3_multigrid.f90.

§ mg3_auto_rstr()

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.

Parameters
[in]box_cChild box
[in,out]box_pParent box
[in]ivIndex of variable
[in]mgMultigrid options

Definition at line 727 of file m_a3_multigrid.f90.

§ mg3_auto_corr()

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.

Parameters
[in,out]box_cChild box
[in]box_pParent box
[in]mgMultigrid options

Definition at line 745 of file m_a3_multigrid.f90.

§ mg3_set_box_tag()

subroutine, public m_a3_multigrid::mg3_set_box_tag ( type(box3_t), intent(inout)  box,
type(mg3_t), intent(in)  mg 
)
Parameters
[in,out]boxBox to operate on
[in]mgMultigrid options

Definition at line 762 of file m_a3_multigrid.f90.

§ mg3_box_corr_lpl()

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 
)
Parameters
[in,out]box_cChild box
[in]box_pParent box
[in]mgMultigrid options

Definition at line 793 of file m_a3_multigrid.f90.

§ mg3_box_gsrb_lpl()

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.

Parameters
[in,out]boxBox to operate on
[in]redblack_cntrIteration counter
[in]mgMultigrid options

Definition at line 803 of file m_a3_multigrid.f90.

§ mg3_box_lpl()

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.

Parameters
[in,out]boxBox to operate on
[in]i_outIndex of variable to store Laplacian in
[in]mgMultigrid options

Definition at line 834 of file m_a3_multigrid.f90.

§ mg3_box_rstr_lpl()

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)

Parameters
[in]box_cChild box to restrict
[in,out]box_pParent box to restrict to
[in]ivVariable to restrict
[in]mgMultigrid options

Definition at line 859 of file m_a3_multigrid.f90.

§ mg3_box_gsrb_lpld()

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.

Parameters
[in,out]boxBox to operate on
[in]redblack_cntrIteration counter
[in]mgMultigrid options

Definition at line 873 of file m_a3_multigrid.f90.

§ mg3_box_lpld()

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.

Parameters
[in,out]boxBox to operate on
[in]i_outIndex of variable to store Laplacian in
[in]mgMultigrid options

Definition at line 910 of file m_a3_multigrid.f90.

§ mg3_box_corr_lpld()

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.

Parameters
[in,out]box_cChild box
[in]box_pParent box
[in]mgMultigrid options

Definition at line 945 of file m_a3_multigrid.f90.

§ mg3_box_corr_lpllsf()

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 
)
Parameters
[in,out]box_cChild box
[in]box_pParent box
[in]mgMultigrid options

Definition at line 1024 of file m_a3_multigrid.f90.

§ mg3_box_gsrb_lpllsf()

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.

Parameters
[in,out]boxBox to operate on
[in]redblack_cntrIteration counter
[in]mgMultigrid options

Definition at line 1077 of file m_a3_multigrid.f90.

§ mg3_box_lpllsf()

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.

Parameters
[in,out]boxBox to operate on
[in]i_outIndex of variable to store Laplacian in
[in]mgMultigrid options

Definition at line 1129 of file m_a3_multigrid.f90.

§ mg3_box_rstr_lpllsf()

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)

Parameters
[in]box_cChild box to restrict
[in,out]box_pParent box to restrict to
[in]ivVariable to restrict
[in]mgMultigrid options

Definition at line 1176 of file m_a3_multigrid.f90.

Variable Documentation

§ mg_normal_box

integer, parameter, public m_a3_multigrid::mg_normal_box = 1

Normal box.

Definition at line 11 of file m_a3_multigrid.f90.

§ mg_lsf_box

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.

§ mg_ceps_box

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.

§ mg_veps_box

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.