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

This module contains the basic types and constants that are used in the 2-dimensional version of Afivo, together with some basic routines. The dimension-independent types and constant are place in m_afivo_types. More...

Data Types

type  a2_cc_methods
 Collection of methods for a cell-centered variable. More...
 
type  a2_loc_t
 Type specifying the location of a cell. More...
 
interface  a2_subr
 Subroutine that gets a box. More...
 
interface  a2_subr_arg
 Subroutine that gets a box and an array of reals. More...
 
interface  a2_subr_bc
 To fill ghost cells near physical boundaries. More...
 
interface  a2_subr_boxes
 Subroutine that gets a list of boxes and a box id. More...
 
interface  a2_subr_boxes_arg
 Subroutine that gets a list of boxes, an id and an array of reals. More...
 
interface  a2_subr_egc
 Subroutine for getting extra ghost cell data (> 1) near physical boundaries. More...
 
interface  a2_subr_prolong
 Subroutine for prolongation. More...
 
interface  a2_subr_rb
 To fill ghost cells near refinement boundaries. More...
 
interface  a2_subr_ref
 Subroutine for setting refinement flags. More...
 
interface  a2_subr_restrict
 Subroutine for restriction. More...
 
type  a2_t
 Type which stores all the boxes and levels, as well as some information about the number of boxes, variables and levels. More...
 
type  box2_t
 The basic building block of afivo: a box with cell-centered and face centered data, and information about its position, neighbors, children etc. More...
 

Functions/Subroutines

subroutine a2_print_info (tree)
 Get tree info. More...
 
pure integer function a2_box_bytes (n_cell, n_var_cell, n_var_face)
 
pure integer function a2_num_boxes_used (tree)
 
pure integer function a2_num_leaves_used (tree)
 
pure integer function a2_num_cells_used (tree)
 
pure real(dp) function a2_total_volume (tree)
 
elemental logical function a2_has_children (box)
 Return .true. if a box has children. More...
 
pure logical function, dimension(size(nbs)) a2_phys_boundary (boxes, id, nbs)
 Return .true. where there is a physical/periodic boundary. Detecting physical boundaries is straightforward (simply test whether neighbors(i) < af_no_box), but periodic boundaries require a comparison of their spatial index. More...
 
pure integer function, dimension(2) a2_get_child_offset (box, nb)
 Get the offset of a box with respect to its parent (e.g. in 2d, there can be a child at offset 0,0, one at n_cell/2,0, one at 0,n_cell/2 and one at n_cell/2, n_cell/2) More...
 
pure integer function, dimension(2) a2_get_ix_on_parent (box, ix)
 Given a cell index on box, get index of the closest cell at its parent. More...
 
pure integer function, dimension(2) a2_get_ix_on_neighb (box, ix, nb)
 Given a cell index on box, get index on a neighbor. More...
 
pure integer function, dimension(2) a2_neighb_offset (nbs)
 Given a list of neighbor directions, compute the index offset. More...
 
pure integer function a2_ix_to_ichild (ix)
 Compute the 'child index' for a box with spatial index ix. With 'child index' we mean the index in the children(:) array of its parent. More...
 
pure integer function, dimension(2) a2_cc_ix (box, r)
 Get the cell index in which r lies. This routine does not check whether r is actually located inside the box. More...
 
pure real(dp) function, dimension(2) a2_r_cc (box, cc_ix)
 Get the location of the cell center with index cc_ix. More...
 
pure real(dp) function, dimension(2) a2_r_fc (box, dim, fc_ix)
 Get the location of the face parallel to dim with index fc_ix. More...
 
pure real(dp) function, dimension(2) a2_rr_cc (box, cc_ix)
 Get a general location with index cc_ix (like a2_r_cc but using reals) More...
 
pure real(dp) function, dimension(2) a2_r_center (box)
 Return the coordinate of the center of a box. More...
 
pure real(dp) function a2_min_dr (tree)
 Return finest dr that is used in the tree. More...
 
pure real(dp) function a2_lvl_dr (tree, lvl)
 Return dr at lvl. More...
 
subroutine a2_set_box_gc (box, nb, iv, gc_scalar, gc_array)
 
pure real(dp) function a2_cyl_radius_cc (box, i)
 Get the radius of the cell center with first index i. More...
 
subroutine a2_cyl_child_weights (box, i, inner, outer)
 Get the normalized weights of the 'inner' and 'outer' children of a cell with index ix. Note that the cell centers of the children are located at -/+ 0.25 dr compared to the parent. More...
 

Variables

integer, parameter a2_num_children = 4
 
integer, parameter a2_child_lowx_lowy = 1
 
integer, parameter a2_child_highx_lowy = 2
 
integer, parameter a2_child_lowx_highy = 3
 
integer, parameter a2_child_highx_highy = 4
 
integer, dimension(2, 4), parameter a2_child_dix = reshape([0,0,1,0,0,1,1,1], [2,4])
 
integer, dimension(2, 4), parameter a2_child_adj_nb = reshape([1,3,2,4,1,2,3,4], [2,4])
 
integer, dimension(2, 4), parameter a2_nb_adj_child = reshape([1,3,2,3,1,4,2,4], [2,4])
 
logical, dimension(2, 4), parameter a2_child_low = reshape([.true., .true., .false., .true., .true., .false., .false., .false.], [2, 4])
 
integer, dimension(2, 4), parameter a2_child_high_01 = reshape([0, 0, 1, 0, 0, 1, 1, 1], [2, 4])
 
integer, parameter a2_num_neighbors = 4
 
integer, parameter a2_neighb_lowx = 1
 
integer, parameter a2_neighb_highx = 2
 
integer, parameter a2_neighb_lowy = 3
 
integer, parameter a2_neighb_highy = 4
 
integer, dimension(2, 4), parameter a2_neighb_dix = reshape([-1,0,1,0,0,-1,0,1], [2,4])
 
logical, dimension(4), parameter a2_neighb_low = [.true., .false., .true., .false.]
 
integer, dimension(2), parameter a2_low_neighbs = [1, 3]
 
integer, dimension(2), parameter a2_high_neighbs = [2, 4]
 
integer, dimension(4), parameter a2_neighb_high_01 = [0, 1, 0, 1]
 
integer, dimension(4), parameter a2_neighb_high_pm = [-1, 1, -1, 1]
 
integer, dimension(4), parameter a2_neighb_rev = [2, 1, 4, 3]
 
integer, dimension(4), parameter a2_neighb_dim = [1, 1, 2, 2]
 

Detailed Description

This module contains the basic types and constants that are used in the 2-dimensional version of Afivo, together with some basic routines. The dimension-independent types and constant are place in m_afivo_types.

Function/Subroutine Documentation

§ a2_print_info()

subroutine m_a2_types::a2_print_info ( type(a2_t), intent(in)  tree)

Get tree info.

Parameters
[in]treeThe tree

Definition at line 221 of file m_a2_types.f90.

§ a2_box_bytes()

pure integer function m_a2_types::a2_box_bytes ( integer, intent(in)  n_cell,
integer, intent(in)  n_var_cell,
integer, intent(in)  n_var_face 
)
Parameters
[in]n_cellnumber of cells per dimension
[in]n_var_cellnumber of cell-centered variables
[in]n_var_facenumber of face-centered variables

Definition at line 252 of file m_a2_types.f90.

§ a2_num_boxes_used()

pure integer function m_a2_types::a2_num_boxes_used ( type(a2_t), intent(in)  tree)

Definition at line 264 of file m_a2_types.f90.

§ a2_num_leaves_used()

pure integer function m_a2_types::a2_num_leaves_used ( type(a2_t), intent(in)  tree)

Definition at line 274 of file m_a2_types.f90.

§ a2_num_cells_used()

pure integer function m_a2_types::a2_num_cells_used ( type(a2_t), intent(in)  tree)

Definition at line 284 of file m_a2_types.f90.

§ a2_total_volume()

pure real(dp) function m_a2_types::a2_total_volume ( type(a2_t), intent(in)  tree)

Definition at line 291 of file m_a2_types.f90.

§ a2_has_children()

elemental logical function m_a2_types::a2_has_children ( type(box2_t), intent(in)  box)

Return .true. if a box has children.

Definition at line 314 of file m_a2_types.f90.

§ a2_phys_boundary()

pure logical function, dimension(size(nbs)) m_a2_types::a2_phys_boundary ( type(box2_t), dimension(:), intent(in)  boxes,
integer, intent(in)  id,
integer, dimension(:), intent(in)  nbs 
)

Return .true. where there is a physical/periodic boundary. Detecting physical boundaries is straightforward (simply test whether neighbors(i) < af_no_box), but periodic boundaries require a comparison of their spatial index.

Parameters
[in]boxesList of boxes
[in]idBox to inspect
[in]nbsNeighbor directions

Definition at line 326 of file m_a2_types.f90.

§ a2_get_child_offset()

pure integer function, dimension(2) m_a2_types::a2_get_child_offset ( type(box2_t), intent(in)  box,
integer, intent(in), optional  nb 
)

Get the offset of a box with respect to its parent (e.g. in 2d, there can be a child at offset 0,0, one at n_cell/2,0, one at 0,n_cell/2 and one at n_cell/2, n_cell/2)

Parameters
[in]boxA child box
[in]nbOptional: get index on parent neighbor

Definition at line 368 of file m_a2_types.f90.

§ a2_get_ix_on_parent()

pure integer function, dimension(2) m_a2_types::a2_get_ix_on_parent ( type(box2_t), intent(in)  box,
integer, dimension(2), intent(in)  ix 
)

Given a cell index on box, get index of the closest cell at its parent.

Parameters
[in]boxA child box
[in]ixIndex on child box

Definition at line 383 of file m_a2_types.f90.

§ a2_get_ix_on_neighb()

pure integer function, dimension(2) m_a2_types::a2_get_ix_on_neighb ( type(box2_t), intent(in)  box,
integer, dimension(2), intent(in)  ix,
integer, intent(in)  nb 
)

Given a cell index on box, get index on a neighbor.

Parameters
[in]boxA box
[in]ixIndex on box
[in]nbNeighbor identifier

Definition at line 391 of file m_a2_types.f90.

§ a2_neighb_offset()

pure integer function, dimension(2) m_a2_types::a2_neighb_offset ( integer, dimension(:), intent(in)  nbs)

Given a list of neighbor directions, compute the index offset.

Parameters
[in]nbsList of neighbor directions

Definition at line 403 of file m_a2_types.f90.

§ a2_ix_to_ichild()

pure integer function m_a2_types::a2_ix_to_ichild ( integer, dimension(2), intent(in)  ix)

Compute the 'child index' for a box with spatial index ix. With 'child index' we mean the index in the children(:) array of its parent.

Parameters
[in]ixSpatial index of the box

Definition at line 416 of file m_a2_types.f90.

§ a2_cc_ix()

pure integer function, dimension(2) m_a2_types::a2_cc_ix ( type(box2_t), intent(in)  box,
real(dp), dimension(2), intent(in)  r 
)

Get the cell index in which r lies. This routine does not check whether r is actually located inside the box.

Definition at line 424 of file m_a2_types.f90.

§ a2_r_cc()

pure real(dp) function, dimension(2) m_a2_types::a2_r_cc ( type(box2_t), intent(in)  box,
integer, dimension(2), intent(in)  cc_ix 
)

Get the location of the cell center with index cc_ix.

Definition at line 432 of file m_a2_types.f90.

§ a2_r_fc()

pure real(dp) function, dimension(2) m_a2_types::a2_r_fc ( type(box2_t), intent(in)  box,
integer, intent(in)  dim,
integer, dimension(2), intent(in)  fc_ix 
)

Get the location of the face parallel to dim with index fc_ix.

Definition at line 440 of file m_a2_types.f90.

§ a2_rr_cc()

pure real(dp) function, dimension(2) m_a2_types::a2_rr_cc ( type(box2_t), intent(in)  box,
real(dp), dimension(2), intent(in)  cc_ix 
)

Get a general location with index cc_ix (like a2_r_cc but using reals)

Definition at line 450 of file m_a2_types.f90.

§ a2_r_center()

pure real(dp) function, dimension(2) m_a2_types::a2_r_center ( type(box2_t), intent(in)  box)

Return the coordinate of the center of a box.

Definition at line 458 of file m_a2_types.f90.

§ a2_min_dr()

pure real(dp) function m_a2_types::a2_min_dr ( type(a2_t), intent(in)  tree)

Return finest dr that is used in the tree.

Returns
Output: dr at the finest lvl of the tree

Definition at line 465 of file m_a2_types.f90.

§ a2_lvl_dr()

pure real(dp) function m_a2_types::a2_lvl_dr ( type(a2_t), intent(in)  tree,
integer, intent(in)  lvl 
)

Return dr at lvl.

Returns
Output: dr at the finest lvl of the tree

Definition at line 472 of file m_a2_types.f90.

§ a2_set_box_gc()

subroutine m_a2_types::a2_set_box_gc ( type(box2_t), intent(inout)  box,
integer, intent(in)  nb,
integer, intent(in)  iv,
real(dp), intent(in), optional  gc_scalar,
real(dp), dimension(box%n_cell), intent(in), optional  gc_array 
)
Parameters
[in,out]boxBox to operate on
[in]nbGhost cell direction
[in]ivGhost cell variable
[in]gc_scalarScalar value for ghost cells
[in]gc_arrayArray for ghost cells

Definition at line 479 of file m_a2_types.f90.

§ a2_cyl_radius_cc()

pure real(dp) function m_a2_types::a2_cyl_radius_cc ( type(box2_t), intent(in)  box,
integer, intent(in)  i 
)

Get the radius of the cell center with first index i.

Definition at line 517 of file m_a2_types.f90.

§ a2_cyl_child_weights()

subroutine m_a2_types::a2_cyl_child_weights ( type(box2_t), intent(in)  box,
integer, intent(in)  i,
real(dp), intent(out)  inner,
real(dp), intent(out)  outer 
)

Get the normalized weights of the 'inner' and 'outer' children of a cell with index ix. Note that the cell centers of the children are located at -/+ 0.25 dr compared to the parent.

Definition at line 527 of file m_a2_types.f90.

Variable Documentation

§ a2_num_children

integer, parameter m_a2_types::a2_num_children = 4

Definition at line 14 of file m_a2_types.f90.

§ a2_child_lowx_lowy

integer, parameter m_a2_types::a2_child_lowx_lowy = 1

Definition at line 15 of file m_a2_types.f90.

§ a2_child_highx_lowy

integer, parameter m_a2_types::a2_child_highx_lowy = 2

Definition at line 16 of file m_a2_types.f90.

§ a2_child_lowx_highy

integer, parameter m_a2_types::a2_child_lowx_highy = 3

Definition at line 17 of file m_a2_types.f90.

§ a2_child_highx_highy

integer, parameter m_a2_types::a2_child_highx_highy = 4

Definition at line 18 of file m_a2_types.f90.

§ a2_child_dix

integer, dimension(2, 4), parameter m_a2_types::a2_child_dix = reshape([0,0,1,0,0,1,1,1], [2,4])

Definition at line 21 of file m_a2_types.f90.

§ a2_child_adj_nb

integer, dimension(2, 4), parameter m_a2_types::a2_child_adj_nb = reshape([1,3,2,4,1,2,3,4], [2,4])

Definition at line 23 of file m_a2_types.f90.

§ a2_nb_adj_child

integer, dimension(2, 4), parameter m_a2_types::a2_nb_adj_child = reshape([1,3,2,3,1,4,2,4], [2,4])

Definition at line 25 of file m_a2_types.f90.

§ a2_child_low

logical, dimension(2, 4), parameter m_a2_types::a2_child_low = reshape([.true., .true., .false., .true., .true., .false., .false., .false.], [2, 4])

Definition at line 27 of file m_a2_types.f90.

§ a2_child_high_01

integer, dimension(2, 4), parameter m_a2_types::a2_child_high_01 = reshape([0, 0, 1, 0, 0, 1, 1, 1], [2, 4])

Definition at line 29 of file m_a2_types.f90.

§ a2_num_neighbors

integer, parameter m_a2_types::a2_num_neighbors = 4

Definition at line 33 of file m_a2_types.f90.

§ a2_neighb_lowx

integer, parameter m_a2_types::a2_neighb_lowx = 1

Definition at line 34 of file m_a2_types.f90.

§ a2_neighb_highx

integer, parameter m_a2_types::a2_neighb_highx = 2

Definition at line 35 of file m_a2_types.f90.

§ a2_neighb_lowy

integer, parameter m_a2_types::a2_neighb_lowy = 3

Definition at line 36 of file m_a2_types.f90.

§ a2_neighb_highy

integer, parameter m_a2_types::a2_neighb_highy = 4

Definition at line 37 of file m_a2_types.f90.

§ a2_neighb_dix

integer, dimension(2, 4), parameter m_a2_types::a2_neighb_dix = reshape([-1,0,1,0,0,-1,0,1], [2,4])

Definition at line 40 of file m_a2_types.f90.

§ a2_neighb_low

logical, dimension(4), parameter m_a2_types::a2_neighb_low = [.true., .false., .true., .false.]

Definition at line 42 of file m_a2_types.f90.

§ a2_low_neighbs

integer, dimension(2), parameter m_a2_types::a2_low_neighbs = [1, 3]

Definition at line 44 of file m_a2_types.f90.

§ a2_high_neighbs

integer, dimension(2), parameter m_a2_types::a2_high_neighbs = [2, 4]

Definition at line 46 of file m_a2_types.f90.

§ a2_neighb_high_01

integer, dimension(4), parameter m_a2_types::a2_neighb_high_01 = [0, 1, 0, 1]

Definition at line 48 of file m_a2_types.f90.

§ a2_neighb_high_pm

integer, dimension(4), parameter m_a2_types::a2_neighb_high_pm = [-1, 1, -1, 1]

Definition at line 50 of file m_a2_types.f90.

§ a2_neighb_rev

integer, dimension(4), parameter m_a2_types::a2_neighb_rev = [2, 1, 4, 3]

Definition at line 53 of file m_a2_types.f90.

§ a2_neighb_dim

integer, dimension(4), parameter m_a2_types::a2_neighb_dim = [1, 1, 2, 2]

Definition at line 55 of file m_a2_types.f90.