Afivo  0.3
Functions/Subroutines
m_af_core Module Reference

This module contains the core routines of Afivo, namely those that deal with initializing and changing the quadtree/octree mesh. More...

Functions/Subroutines

subroutine, public af_add_cc_variable (tree, name, write_out, n_copies, ix, write_binary)
 Add cell-centered variable. More...
 
subroutine, public af_add_fc_variable (tree, name, ix, write_binary)
 Add face-centered variable. More...
 
integer function, public af_find_cc_variable (tree, name)
 Find index of cell-centered variable. More...
 
integer function, public af_find_fc_variable (tree, name)
 Find index of face-centered variable. More...
 
subroutine, public af_init (tree, n_cell, r_max, grid_size, periodic, r_min, coord, mem_limit_gb)
 Initialize a NDIM-d octree/quadtree grid. More...
 
subroutine af_set_coarse_grid (tree, coarse_grid_size, periodic_dims)
 Create the coarse grid. More...
 
subroutine, public af_set_cc_methods (tree, iv, bc, rb, prolong, restrict, bc_custom, funcval, prolong_limiter)
 Set the methods for a cell-centered variable. More...
 
subroutine, public af_destroy (tree)
 "Destroy" the data in a tree. Since we don't use pointers, you can also just let a tree get out of scope More...
 
subroutine create_index_array (nx, periodic, id_array)
 Create an array for easy lookup of indices. More...
 
subroutine set_leaves_parents (boxes, level)
 Create a list of leaves and a list of parents for a level. More...
 
subroutine, public af_init_box (tree, id)
 Mark box as active and allocate data storage for a box, for its cell- and face-centered data. More...
 
subroutine af_deactivate_box (box)
 Mark box as inactive, but keep storage for cell- and face-centered data to avoid reallocating this. More...
 
subroutine set_neighbs (boxes, id)
 
integer function find_neighb (boxes, id, dix)
 Get the id of all neighbors of boxes(id), through its parent. More...
 
subroutine, public af_refine_up_to_lvl (tree, lvl)
 Refine a tree up to a given refinement lvl. More...
 
subroutine, public af_adjust_refinement (tree, ref_subr, ref_info, ref_buffer, ref_links)
 Adjust the refinement of a tree using the user-supplied ref_subr. The optional argument ref_buffer controls over how many cells neighbors are affected by refinement flags. More...
 
subroutine auto_restrict (tree, id)
 Try to automatically restrict to box with index id. More...
 
subroutine auto_prolong (tree, ref_info)
 Try to automatically prolong to all new boxes. More...
 
subroutine get_free_ids (tree, ids)
 Get free ids from the boxes(:) array to store new boxes in. These ids are always consecutive. More...
 
subroutine consistent_ref_flags (tree, ref_flags, ref_subr, ref_buffer, ref_links)
 Given the refinement function, return consistent refinement flags, that ensure that the tree is still balanced. Furthermore, it cannot derefine the base level, and it cannot refine above af_max_lvl. The argument ref_flags is changed: for boxes that will be refined it holds af_refine, for boxes that will be derefined it holds af_derefine. More...
 
subroutine ensure_two_one_balance (tree, ref_flags)
 Adjust refinement flags to ensure 2-1 balance is maintained. More...
 
subroutine handle_derefinement_flags (tree, ref_flags)
 
subroutine cell_to_ref_flags (cell_flags, nc, ref_flags, tree, id, ref_buffer)
 Given the cell refinement flags of a box, set the refinement flag for that box and potentially also its neighbors (in case of refinement near a boundary). More...
 
subroutine remove_children (tree, id)
 Remove the children of box id. More...
 
subroutine add_children (tree, id, c_ids)
 Add children to box id, using the indices in c_ids. More...
 
subroutine set_child_ids (p_ids, c_ids, boxes)
 Create a list c_ids(:) of all the children of p_ids(:). This is used after a level has been refined. More...
 
subroutine, public af_consistent_fluxes (tree, f_ixs)
 Restrict fluxes from children to parents on refinement boundaries. More...
 
subroutine flux_from_children (boxes, id, nb, f_ixs)
 The neighbor nb has no children and id does, so set flux on the neighbor from our children. This ensures flux consistency at refinement boundary. More...
 

Detailed Description

This module contains the core routines of Afivo, namely those that deal with initializing and changing the quadtree/octree mesh.

Function/Subroutine Documentation

◆ af_add_cc_variable()

subroutine, public m_af_core::af_add_cc_variable ( type(af_t), intent(inout)  tree,
character(len=*), intent(in)  name,
logical, intent(in), optional  write_out,
integer, intent(in), optional  n_copies,
integer, intent(out), optional  ix,
logical, intent(in), optional  write_binary 
)

Add cell-centered variable.

Todo:
ix as third argument?
Parameters
[in,out]treeTree to add variable to
[in]nameName of the variable
[in]write_outInclude variable in output
[in]write_binaryInclude variable in binary output (for restarting)
[in]n_copiesHow many copies of variable to store (default: 1)
[out]ixOn output: index of variable

Definition at line 26 of file m_af_core.f90.

◆ af_add_fc_variable()

subroutine, public m_af_core::af_add_fc_variable ( type(af_t), intent(inout)  tree,
character(len=*), intent(in)  name,
integer, intent(out), optional  ix,
logical, intent(in), optional  write_binary 
)

Add face-centered variable.

Parameters
[in,out]treeTree to add variable to
[in]nameName of the variable
[out]ixOn output: index of variable
[in]write_binaryInclude variable in binary output

Definition at line 76 of file m_af_core.f90.

◆ af_find_cc_variable()

integer function, public m_af_core::af_find_cc_variable ( type(af_t), intent(in)  tree,
character(len=*), intent(in)  name 
)

Find index of cell-centered variable.

Definition at line 102 of file m_af_core.f90.

◆ af_find_fc_variable()

integer function, public m_af_core::af_find_fc_variable ( type(af_t), intent(in)  tree,
character(len=*), intent(in)  name 
)

Find index of face-centered variable.

Definition at line 120 of file m_af_core.f90.

◆ af_init()

subroutine, public m_af_core::af_init ( type(af_t), intent(inout)  tree,
integer, intent(in)  n_cell,
real(dp), dimension(ndim), intent(in)  r_max,
integer, dimension(ndim), intent(in)  grid_size,
logical, dimension(ndim), intent(in), optional  periodic,
real(dp), dimension(ndim), intent(in), optional  r_min,
integer, intent(in), optional  coord,
real(dp), intent(in), optional  mem_limit_gb 
)

Initialize a NDIM-d octree/quadtree grid.

Parameters
[in,out]treeThe tree to initialize
[in]n_cellBoxes have n_cell^dim cells
[in]r_maxMaximal coordinates of the domain
[in]grid_sizeSize of the coarse grid
[in]periodicTrue for periodic dimensions
[in]r_minLowest coordinate, default is (0., 0., 0.)
[in]coordSelect coordinate type
[in]mem_limit_gbMemory limit in GByte

Definition at line 138 of file m_af_core.f90.

◆ af_set_coarse_grid()

subroutine m_af_core::af_set_coarse_grid ( type(af_t), intent(inout)  tree,
integer, dimension(ndim), intent(in)  coarse_grid_size,
logical, dimension(ndim), intent(in), optional  periodic_dims 
)
private

Create the coarse grid.

Parameters
[in,out]treeTree for which we set the base
[in]coarse_grid_sizeSize of coarse grid (in cells)
[in]periodic_dimsWhether dimensions are periodic (default: false)

Definition at line 198 of file m_af_core.f90.

◆ af_set_cc_methods()

subroutine, public m_af_core::af_set_cc_methods ( type(af_t), intent(inout)  tree,
integer, intent(in)  iv,
procedure(af_subr_bc), optional  bc,
procedure(af_subr_rb), optional  rb,
procedure(af_subr_prolong), optional  prolong,
procedure(af_subr_restrict), optional  restrict,
procedure(af_subr_bc_custom), optional  bc_custom,
procedure(af_subr_funcval), optional  funcval,
integer, intent(in), optional  prolong_limiter 
)

Set the methods for a cell-centered variable.

Parameters
[in,out]treeTree to operate on
[in]ivIndex of variable
bcBoundary condition method
rbRefinement boundary method
prolongProlongation method
restrictRestriction method
bc_customCustom b.c. method
funcvalVariable defined by function Type of limiter to use for prolongation (of values or ghost cells)

Definition at line 335 of file m_af_core.f90.

◆ af_destroy()

subroutine, public m_af_core::af_destroy ( type(af_t), intent(out)  tree)

"Destroy" the data in a tree. Since we don't use pointers, you can also just let a tree get out of scope

Definition at line 423 of file m_af_core.f90.

◆ create_index_array()

subroutine m_af_core::create_index_array ( integer, dimension(ndim), intent(in)  nx,
logical, dimension(ndim), intent(in)  periodic,
integer, dimension(dtimes(:)), intent(inout), allocatable  id_array 
)
private

Create an array for easy lookup of indices.

Definition at line 428 of file m_af_core.f90.

◆ set_leaves_parents()

subroutine m_af_core::set_leaves_parents ( type(box_t), dimension(:), intent(in)  boxes,
type(lvl_t), intent(inout)  level 
)
private

Create a list of leaves and a list of parents for a level.

Parameters
[in]boxesList of boxes
[in,out]levelLevel type which contains the indices of boxes

Definition at line 496 of file m_af_core.f90.

◆ af_init_box()

subroutine, public m_af_core::af_init_box ( type(af_t), intent(inout)  tree,
integer, intent(in)  id 
)

Mark box as active and allocate data storage for a box, for its cell- and face-centered data.

Parameters
[in]idBox id

Definition at line 531 of file m_af_core.f90.

◆ af_deactivate_box()

subroutine m_af_core::af_deactivate_box ( type(box_t), intent(inout)  box)
private

Mark box as inactive, but keep storage for cell- and face-centered data to avoid reallocating this.

Definition at line 575 of file m_af_core.f90.

◆ set_neighbs()

subroutine m_af_core::set_neighbs ( type(box_t), dimension(:), intent(inout)  boxes,
integer, intent(in)  id 
)
private

Definition at line 587 of file m_af_core.f90.

◆ find_neighb()

integer function m_af_core::find_neighb ( type(box_t), dimension(:), intent(in)  boxes,
integer, intent(in)  id,
integer, dimension(ndim), intent(in)  dix 
)
private

Get the id of all neighbors of boxes(id), through its parent.

Parameters
[in]boxesList with all the boxes
[in]idBox whose neighbor we are looking for

Definition at line 628 of file m_af_core.f90.

◆ af_refine_up_to_lvl()

subroutine, public m_af_core::af_refine_up_to_lvl ( type(af_t), intent(inout)  tree,
integer, intent(in)  lvl 
)

Refine a tree up to a given refinement lvl.

Parameters
[in,out]treeThe tree to adjust
[in]lvlRefine up to this lvl

Definition at line 656 of file m_af_core.f90.

◆ af_adjust_refinement()

subroutine, public m_af_core::af_adjust_refinement ( type(af_t), intent(inout)  tree,
procedure(af_subr_ref ref_subr,
type(ref_info_t), intent(inout)  ref_info,
integer, intent(in), optional  ref_buffer,
integer, dimension(:, :), intent(in), optional  ref_links 
)

Adjust the refinement of a tree using the user-supplied ref_subr. The optional argument ref_buffer controls over how many cells neighbors are affected by refinement flags.

On input, the tree should be balanced. On output, the tree is still balanced, and its refinement is updated (with at most one level per call).

Parameters
[in,out]treeThe tree to adjust
ref_subrRefinement function
[in,out]ref_infoInformation about refinement
[in]ref_bufferBuffer width (in cells)
[in]ref_linksLists of linked boxes which should have the same refinement

Definition at line 689 of file m_af_core.f90.

◆ auto_restrict()

subroutine m_af_core::auto_restrict ( type(af_t), intent(inout)  tree,
integer, intent(in)  id 
)
private

Try to automatically restrict to box with index id.

Definition at line 817 of file m_af_core.f90.

◆ auto_prolong()

subroutine m_af_core::auto_prolong ( type(af_t), intent(inout)  tree,
type(ref_info_t), intent(in)  ref_info 
)
private

Try to automatically prolong to all new boxes.

Definition at line 835 of file m_af_core.f90.

◆ get_free_ids()

subroutine m_af_core::get_free_ids ( type(af_t), intent(inout)  tree,
integer, dimension(:), intent(out)  ids 
)

Get free ids from the boxes(:) array to store new boxes in. These ids are always consecutive.

Parameters
[out]idsArray which will be filled with free box ids
Todo:
when doing AMR in parallel, perhaps move some of the code outside the critical construct

Definition at line 877 of file m_af_core.f90.

◆ consistent_ref_flags()

subroutine m_af_core::consistent_ref_flags ( type(af_t), intent(inout)  tree,
integer, dimension(:), intent(inout)  ref_flags,
procedure(af_subr_ref ref_subr,
integer, intent(in)  ref_buffer,
integer, dimension(:, :), intent(in), optional  ref_links 
)
private

Given the refinement function, return consistent refinement flags, that ensure that the tree is still balanced. Furthermore, it cannot derefine the base level, and it cannot refine above af_max_lvl. The argument ref_flags is changed: for boxes that will be refined it holds af_refine, for boxes that will be derefined it holds af_derefine.

Parameters
[in,out]treeTree for which we set refinement flags
[in,out]ref_flagsList of refinement flags for all boxes(:)
ref_subrUser-supplied refinement function.
[in]ref_bufferBuffer width (in cells)
[in]ref_linksLists of linked boxes which should have the same refinement

Definition at line 921 of file m_af_core.f90.

◆ ensure_two_one_balance()

subroutine m_af_core::ensure_two_one_balance ( type(af_t), intent(inout)  tree,
integer, dimension(:), intent(inout)  ref_flags 
)

Adjust refinement flags to ensure 2-1 balance is maintained.

Definition at line 1008 of file m_af_core.f90.

◆ handle_derefinement_flags()

subroutine m_af_core::handle_derefinement_flags ( type(af_t), intent(inout)  tree,
integer, dimension(:), intent(inout)  ref_flags 
)
private

Definition at line 1051 of file m_af_core.f90.

◆ cell_to_ref_flags()

subroutine m_af_core::cell_to_ref_flags ( integer, dimension(dtimes(nc)), intent(in)  cell_flags,
integer, intent(in)  nc,
integer, dimension(:), intent(inout)  ref_flags,
type(af_t), intent(in)  tree,
integer, intent(in)  id,
integer, intent(in)  ref_buffer 
)
private

Given the cell refinement flags of a box, set the refinement flag for that box and potentially also its neighbors (in case of refinement near a boundary).

Parameters
[in]ncn_cell for the box
[in]cell_flagsCell refinement flags
[in,out]ref_flagsBox refinement flags for this thread
[in]treeFull tree
[in]idWhich box is considered
[in]ref_bufferBuffer cells around refinement

Definition at line 1087 of file m_af_core.f90.

◆ remove_children()

subroutine m_af_core::remove_children ( type(af_t), intent(inout)  tree,
integer, intent(in)  id 
)

Remove the children of box id.

Parameters
[in]idId of box whose children will be removed

Definition at line 1143 of file m_af_core.f90.

◆ add_children()

subroutine m_af_core::add_children ( type(af_t), intent(inout)  tree,
integer, intent(in)  id,
integer, dimension(af_num_children), intent(in)  c_ids 
)
private

Add children to box id, using the indices in c_ids.

Parameters
[in]idId of box that gets children
[in]c_idsFree ids for the children

Definition at line 1180 of file m_af_core.f90.

◆ set_child_ids()

subroutine m_af_core::set_child_ids ( integer, dimension(:), intent(in)  p_ids,
integer, dimension(:), intent(inout), allocatable  c_ids,
type(box_t), dimension(:), intent(in)  boxes 
)
private

Create a list c_ids(:) of all the children of p_ids(:). This is used after a level has been refined.

Parameters
[in]p_idsAll the parents ids
[in,out]c_idsOutput: all the children's ids
[in]boxesList of all the boxes

Definition at line 1229 of file m_af_core.f90.

◆ af_consistent_fluxes()

subroutine, public m_af_core::af_consistent_fluxes ( type(af_t), intent(inout)  tree,
integer, dimension(:), intent(in)  f_ixs 
)

Restrict fluxes from children to parents on refinement boundaries.

Parameters
[in,out]treeTree to operate on
[in]f_ixsIndices of the fluxes

Definition at line 1249 of file m_af_core.f90.

◆ flux_from_children()

subroutine m_af_core::flux_from_children ( type(box_t), dimension(:), intent(inout)  boxes,
integer, intent(in)  id,
integer, intent(in)  nb,
integer, dimension(:), intent(in)  f_ixs 
)
private

The neighbor nb has no children and id does, so set flux on the neighbor from our children. This ensures flux consistency at refinement boundary.

Parameters
[in,out]boxesList of all the boxes
[in]idId of box for which we set fluxes
[in]nbDirection in which fluxes are set
[in]f_ixsIndices of the fluxes

Definition at line 1278 of file m_af_core.f90.