FVBCs System

For an overview of MOOSE FV please see Finite Volume Design Decisions in MOOSE.

The finite volume method (FVM) distinguishes between two types of boundary conditions.

* FVDirichletBC prescribes values of the FVM variables on the boundary. This boundary condition acts similarly to Dirichlet boundary conditions in FEM but it is implemented using a ghost element method.

* FVFluxBC prescribes the flux on a boundary. This boundary condition is similar to integrated boundary conditions in FEM.

Currently, the FVDirichletBC category only contains a single class that applies a fixed value on the boundary. In the future, more specialized classes will be added.

FVBCs block

FVM boundary conditions are added to simulation input files in the FVBCs as in the example below.

Listing 1: Example of the FVBCs block in a MOOSE input file.

[FVBCs]
  [left]
    type = FVNeumannBC
    variable = v
    boundary = left
    value = 5
  []
  [right]
    type = FVDirichletBC
    variable = v
    boundary = right
    value = 42
  []
[]
(moose/test/tests/fvkernels/fv_simple_diffusion/neumann.i)

In this example input, a diffusion equation with flux boundary conditions on the left and Dirichlet boundary conditions on the right is solved. To understand the differences between these two boundary conditions, let's start with the diffusion equation:

and the boundary conditions on the left:

where is the outward normal and on the right:

For seeing how the flux boundary condition is applied, the diffusion equation is integrated over the extent of an element adjacent to the left boundary and Gauss' theorem is applied to the divergence:

where is the element volume, are all faces that belong to the left sideset, are all faces, and is the area of face. Flux boundary conditions are applied by replacing appropriate terms in the FVM balance by the value of the flux prescribed on the boundary.

Dirichlet boundary conditions are applied differently. Let us first write a balance equation for an element that is adjacent to the right boundary:

MOOSE uses the ghost element method to apply Dirichlet boundary conditions for FVM. The process of using a ghost elements is the following:

  1. Place a virtual element across the Dirichlet boundary.

  2. Compute the value of in the ghost element as the extrapolation of the element value and the desired value on the boundary.

  3. Assign the value of in the ghost element.

  4. Evaluate the numerical fluxes as if you were on an interior face.

For implementing the ghost element method an extrapolation must be selected. Currently, MOOSE FVM only supports linear extrapolation. If the value of the Dirichlet boundary condition is denoted by and the value in the element is denoted by , then the ghost element value is:

The parameters available in boundary conditions are equivalent to FEM boundary conditions and are not discussed in detail here.

FVBCs source code: FVDirichletBC

FVDirichletBC objects assigns a constant value on a boundary. Implementation of a FVM Dirichlet boundary condition usually only requires overriding the boundaryValue method. The boundaryValue method must return the value of the variable on the Dirichlet boundary.

Listing 2: Example source code for FVDirichletBC.

#include "FVDirichletBC.h"

registerMooseObject("MooseApp", FVDirichletBC);

InputParameters
FVDirichletBC::validParams()
{
  InputParameters params = FVDirichletBCBase::validParams();
  params.addClassDescription("Defines a Dirichlet boundary condition for finite volume method.");
  params.addRequiredParam<Real>("value", "value to enforce at the boundary face");
  return params;
}

FVDirichletBC::FVDirichletBC(const InputParameters & parameters)
  : FVDirichletBCBase(parameters), _val(getParam<Real>("value"))
{
}

ADReal
FVDirichletBC::boundaryValue(const FaceInfo & /*fi*/, const Moose::StateArg & /*state*/) const
{
  return _val;
}
(moose/framework/src/fvbcs/FVDirichletBC.C)

FVBCs source code: FVFluxBC

FVNeumannBC objects assign a constant flux on a boundary. Implementation of a flux boundary condition usually only requires overriding the computeQpResidual() method. FVNeumannBC reads a constant value from the parameters and then returns it in computeQpResidual().

Listing 3: Example source code for FVNeumannBC.

#include "FVNeumannBC.h"

registerMooseObject("MooseApp", FVNeumannBC);

InputParameters
FVNeumannBC::validParams()
{
  InputParameters params = FVFluxBC::validParams();
  params.addClassDescription("Neumann boundary condition for finite volume method.");
  params.addParam<Real>("value", 0.0, "The value of the flux crossing the boundary.");
  return params;
}

FVNeumannBC::FVNeumannBC(const InputParameters & parameters)
  : FVFluxBC(parameters), _value(getParam<Real>("value"))
{
}

ADReal
FVNeumannBC::computeQpResidual()
{
  return -_value;
}
(moose/framework/src/fvbcs/FVNeumannBC.C)

FVBCs source code: FVBurgersOutflowBC

Flux boundary conditions can be more complicated than assigning a constant value. In this example, the outflow contribution for the Burgers' equation is computed. The relevant term is (note 1D):

where and are the values of on the left and right boundaries of the element and and are the outward normals on these faces. Let's assume that the left side is a boundary face where the FVBurgersOutflowBC is applied. On that boundary we have . The FVBurgersOutflowBC boundary condition uses upwinding, i.e. it uses the element value as boundary values on outflow faces.

The code listed below first checks if the left is actually an outflow side by using the cell value of the (element average, upwinding!) and dotting it with the normal. If , then the left is an outflow side. In this case the contribution is added, otherwise no contribution is added.

Listing 4: Outflow boundary condition for the Burgers' equation.

FVBurgersOutflowBC::computeQpResidual()
{
  mooseAssert(_face_info->elem().dim() == 1, "FVBurgersOutflowBC works only in 1D");

  ADReal r = 0;
  // only add this on outflow faces
  if (_u[_qp] * _normal(0) > 0)
    r = 0.5 * _u[_qp] * _u[_qp] * _normal(0);
  return r;
}
(moose/test/src/fvbcs/FVBurgersOutflowBC.C)

In this case, the boundary condition does not represent a fixed inflow, but rather a free outflow condition.

Available Objects

  • Moose App
  • FVADFunctorDirichletBCUses the value of a functor to set a Dirichlet boundary value.
  • FVBoundaryIntegralValueConstraintThis class is used to enforce integral of phi = boundary area * phi_0 with a Lagrange multiplier approach.
  • FVConstantScalarOutflowBCConstant velocity scalar advection boundary conditions for finite volume method.
  • FVDirichletBCDefines a Dirichlet boundary condition for finite volume method.
  • FVFunctionDirichletBCImposes the essential boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
  • FVFunctionNeumannBCNeumann boundary condition for finite volume method.
  • FVFunctorDirichletBCUses the value of a functor to set a Dirichlet boundary value.
  • FVFunctorNeumannBCNeumann boundary condition for the finite volume method.
  • FVNeumannBCNeumann boundary condition for finite volume method.
  • FVOrthogonalBoundaryDiffusionImposes an orthogonal diffusion boundary term with specified boundary function.
  • FVPostprocessorDirichletBCDefines a Dirichlet boundary condition for finite volume method.