Afivo  0.3
Boundary conditions

Introduction

Unless a computational domain is fully periodic, there will be physical boundaries. Afivo implements boundary conditions using ghost cells, see Filling ghost cells. Currently, Dirichlet and Neumann type boundary conditions are supported. An example of how to fill boundary conditions is given in boundary_conditions.f90.

Writing a boundary condition routine

The user has to provide routines that will be called near physical boundaries when adding variables through m_af_core::af_add_cc_variable(). The interface of these routines is given by m_af_types::af_subr_bc. An example is shown below:

subroutine boundary_method(box, nb, iv, coords, bc_val, bc_type)
type(box_t), intent(in) :: box
integer, intent(in) :: nb
integer, intent(in) :: iv
real(dp), intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
real(dp), intent(out) :: bc_val(box%n_cell**(NDIM-1))
integer, intent(out) :: bc_type
integer :: nc
nc = box%n_cell
! Below the solution is specified in the approriate ghost cells
select case (nb)
case (af_neighb_lowx) ! Lower-x direction
bc_type = af_bc_dirichlet
bc_val = 1.0_dp
case default
bc_type = af_bc_neumann
bc_val = 0.0_dp
end select
end subroutine boundary_method