FVKernels System
For an overview of MOOSE FV please see Finite Volume Design Decisions in MOOSE.
For the finite volume method (FVM), FVKernels
are the base class for FVFluxKernel
, FVElementalKernel
. These specialized objects satisfy the following tasks:
* FVFluxKernel
represents numerical fluxes evaluate on the element faces. These terms originate from applying Gauss' divergence theorem.
* FVElementalKernel
represents terms that do not contain a spatial derivative so that Gauss' theorem cannot be applied. These terms include time derivatives, externally imposed source terms, and reaction terms.
Note: Currently, the FVElementalKernel
category only contains kernels (subclasses) representing time derivatives. Kernels representing externally imposed sources or reaction terms will be added in the near future.
In the documentation that follows, we will use '-' and '' to represent different sides of a face. This is purely notation. In the MOOSE code base, the '-' side is represented with an _elem
suffix and the '' side is represented with a _neighbor
suffix. We could just as well have chosen _left
and _right
, or _1
and _2
, or _minus
and _plus
, but for consistency with previous MOOSE framework code such as discontinuous Galerkin kernels and node-face constraints, we have elected to go with the _elem
and _neighbor
suffixes.
FVKernels block
FVM kernels are added to simulation input files in the FVKernels
block. The FVKernels
block in the example below sets up a transient diffusion problem defined by the equation:
The time derivative term corresponds to the kernel named time
, while the diffusion term is represented by the kernel named diff
.
[FVKernels]
[./time]
type = FVFunctorTimeKernel
variable = v
[../]
[diff]
type = FVDiffusion
variable = v
coeff = coeff
[]
[]
(moose/test/tests/fvkernels/fv_simple_diffusion/transient.i)The FVTimeKernel
in the example derives from FVElementalKernel
so it's a volumetric contribution to the residual, while the FVDiffusion
kernel is an FVFluxKernel
and it's a face contribution to the residual. The remaining MOOSE syntax is what you would expect to see in finite element kernel objects:
* variable
refers to the variable that this kernel is acting on (i.e. into which equation does the residual of this term go). This must be a finite-volume variable (defined with fv = true
) for all FVM kernels.
* coeff
in kernel diff
is a material property corresponding to the heat conduction or diffusion coefficient.
The next example shows an FVKernels
block that solves the one-dimensional Burgers' equation. The Burgers' equation for speed v
is given by:
[FVKernels]
[./burgers]
type = FVBurgers1D
variable = v
[../]
[./time]
type = FVTimeKernel
variable = v
[../]
[]
(moose/test/tests/fvkernels/fv_burgers/fv_burgers.i)Note that the FVBurgers1D
kernel only works for one-dimensional problems. In this example, the exact same time derivative kernels as for the diffusion example is used, but the spatial derivative term is different.
Boundary conditions are not discussed in these examples. Look at syntax files for details about boundary conditions.
Available Objects
- Moose App
- FVAdvectionResidual contribution from advection operator for finite volume method.
- FVAnisotropicDiffusionComputes residual for anisotropic diffusion operator for finite volume method.
- FVBodyForceDemonstrates the multiple ways that scalar values can be introduced into finite volume kernels, e.g. (controllable) constants, functions, and postprocessors.
- FVBoundedValueConstraintThis class is used to enforce a min or max value for a finite volume variable
- FVCoupledForceImplements a source term proportional to the value of a coupled variable.
- FVDiffusionComputes residual for diffusion operator for finite volume method.
- FVDivergenceComputes the residual coming from the divergence of a vector fieldthat can be represented as a functor.
- FVFunctorTimeKernelResidual contribution from time derivative of an AD functor (default is the variable this kernel is acting upon if the 'functor' parameter is not supplied) for the finite volume method.
- FVIntegralValueConstraintThis class is used to enforce integral of phi = volume * phi_0 with a Lagrange multiplier approach.
- FVMassMatrixComputes a 'mass matrix', which will just be a diagonal matrix for the finite volume method, meant for use in preconditioning schemes which require one
- FVMatAdvectionComputes the residual of advective term using finite volume method.
- FVOrthogonalDiffusionImposes an orthogonal diffusion term.
- FVPointValueConstraintThis class is used to enforce integral of phi = volume * phi_0 with a Lagrange multiplier approach.
- FVReactionSimple consuming reaction term
- FVScalarLagrangeMultiplierThis class is used to enforce integral of phi = volume * phi_0 with a Lagrange multiplier approach.
- FVTimeKernelResidual contribution from time derivative of a variable for the finite volume method.
FVKernel source code: FVDiffusion example
First, FVFluxKernels
are discussed. FVFluxKernels
are used to calculate numerical flux contributions from face (surface integral) terms to the residual. The residual contribution is implemented by overriding the computeQpResidual
function.
In the FVM, one solves for the averages of the variables over each element. The values of the variables on the faces are unknown and must be computed from the cell average values. This interpolation/reconstruction determines the accuracy of the FVM. The discussion is based on the example of FVDiffusion
that discretizes the diffusion term using a central difference approximation.
#include "FVDiffusion.h"
registerMooseObject("MooseApp", FVDiffusion);
InputParameters
FVDiffusion::validParams()
{
InputParameters params = FVFluxKernel::validParams();
params.addClassDescription("Computes residual for diffusion operator for finite volume method.");
params.addRequiredParam<MooseFunctorName>("coeff", "diffusion coefficient");
MooseEnum coeff_interp_method("average harmonic", "harmonic");
params.addParam<MooseEnum>(
"coeff_interp_method",
coeff_interp_method,
"Switch that can select face interpolation method for diffusion coefficients.");
params.set<unsigned short>("ghost_layers") = 2;
return params;
}
FVDiffusion::FVDiffusion(const InputParameters & params)
: FVFluxKernel(params),
_coeff(getFunctor<ADReal>("coeff")),
_coeff_interp_method(
Moose::FV::selectInterpolationMethod(getParam<MooseEnum>("coeff_interp_method")))
{
if ((_var.faceInterpolationMethod() == Moose::FV::InterpMethod::SkewCorrectedAverage) &&
(_tid == 0))
adjustRMGhostLayers(std::max((unsigned short)(3), _pars.get<unsigned short>("ghost_layers")));
}
ADReal
FVDiffusion::computeQpResidual()
{
using namespace Moose::FV;
const auto state = determineState();
auto dudn = gradUDotNormal(state);
ADReal coeff;
// If we are on internal faces, we interpolate the diffusivity as usual
if (_var.isInternalFace(*_face_info))
{
const ADReal coeff_elem = _coeff(elemArg(), state);
const ADReal coeff_neighbor = _coeff(neighborArg(), state);
// If the diffusion coefficients are zero, then we can early return 0 (and avoid warnings if we
// have a harmonic interpolation)
if (!coeff_elem.value() && !coeff_neighbor.value())
return 0;
interpolate(_coeff_interp_method, coeff, coeff_elem, coeff_neighbor, *_face_info, true);
}
// Else we just use the boundary values (which depend on how the diffusion
// coefficient is constructed)
else
{
const auto face = singleSidedFaceArg();
coeff = _coeff(face, state);
}
return -1 * coeff * dudn;
}
(moose/framework/src/fvkernels/FVDiffusion.C)The kernel FVDiffusion
discretizes the diffusion term . Integrating over the extend of an element and using Gauss' theorem leads to:
The term in parenthesis in the surface integral on the right hand side must be implemented in the FVKernel
. However, there is one more step before we can implement the kernel. We must determine how the values of and depend on the values of and on the '+' and '-' side of the face . In this example, the following approximation is used:
This is a central difference approximation of the gradient on the face that neglects cross diffusion terms.
Now, the implementation of this numerical flux into FVDiffusion::computeQpResidual
is discussed.
* the kernel provides the '-' and '+' values of the variable as _u_elem[_qp]
and _u_neighbor[_qp]
* the values of the material properties on the '-' side of the face is obtained by _coeff_elem(getADMaterialProperty<Real>("coeff"))
while the '+' side value is obtained by calling getNeighborADMaterialProperty<Real>("coeff")
.
* geometric information about the '-' and '+' adjacent elements is available from the face_info
object.
The implementation is then straight forward. The first line of the code computes dudn
which corresponds to the term:
while the second line computes k
corresponding to:
The minus sign originates from the minus sign in the original expression. Flow from '-' to '+ is defined to be positive.
FVKernel source code: FVMatAdvection example
In this example the advection term:
is discretized using upwinding. The velocity is denoted by and represents a passive scalar quantity advected by the flow. Upwinding is a strategy that approximates the value of a variable on a face by taking the value from the upwind element (i.e. the element where the flow originates from).
#include "FVDiffusion.h"
registerMooseObject("MooseApp", FVDiffusion);
InputParameters
FVDiffusion::validParams()
{
InputParameters params = FVFluxKernel::validParams();
params.addClassDescription("Computes residual for diffusion operator for finite volume method.");
params.addRequiredParam<MooseFunctorName>("coeff", "diffusion coefficient");
MooseEnum coeff_interp_method("average harmonic", "harmonic");
params.addParam<MooseEnum>(
"coeff_interp_method",
coeff_interp_method,
"Switch that can select face interpolation method for diffusion coefficients.");
params.set<unsigned short>("ghost_layers") = 2;
return params;
}
FVDiffusion::FVDiffusion(const InputParameters & params)
: FVFluxKernel(params),
_coeff(getFunctor<ADReal>("coeff")),
_coeff_interp_method(
Moose::FV::selectInterpolationMethod(getParam<MooseEnum>("coeff_interp_method")))
{
if ((_var.faceInterpolationMethod() == Moose::FV::InterpMethod::SkewCorrectedAverage) &&
(_tid == 0))
adjustRMGhostLayers(std::max((unsigned short)(3), _pars.get<unsigned short>("ghost_layers")));
}
ADReal
FVDiffusion::computeQpResidual()
{
using namespace Moose::FV;
const auto state = determineState();
auto dudn = gradUDotNormal(state);
ADReal coeff;
// If we are on internal faces, we interpolate the diffusivity as usual
if (_var.isInternalFace(*_face_info))
{
const ADReal coeff_elem = _coeff(elemArg(), state);
const ADReal coeff_neighbor = _coeff(neighborArg(), state);
// If the diffusion coefficients are zero, then we can early return 0 (and avoid warnings if we
// have a harmonic interpolation)
if (!coeff_elem.value() && !coeff_neighbor.value())
return 0;
interpolate(_coeff_interp_method, coeff, coeff_elem, coeff_neighbor, *_face_info, true);
}
// Else we just use the boundary values (which depend on how the diffusion
// coefficient is constructed)
else
{
const auto face = singleSidedFaceArg();
coeff = _coeff(face, state);
}
return -1 * coeff * dudn;
}
(moose/framework/src/fvkernels/FVDiffusion.C)Integrating the advection term over the element and using Gauss' theorem leads to:
This term in parenthesis on the right hand side is approximated using upwinding:
where
and if and otherwise. By convention, the normal points from the '-' side to the '+' side.
The implementation is straight forward. In the constructor the '-' and '' velocities are obtained as RealVectorValue
material properties. The average is computed and stored in variable v_avg
. The direction of the flow is determined using the inner product of v_avg * _normal
and the residual is then computed using either the '-' value of given by _u_elem[_qp]
or the '' value given by _u_neighbor [_qp]
.
FVKernel source code: FVTimeKernel
This example demonstrates source code for an FVElementalKernel
. FVElementalKernel
are volumetric terms. In this case, the kernel is FVTimeKernel
.
FVTimeKernel::computeQpResidual()
{
return _u_dot[_qp];
}
(moose/framework/src/fvkernels/FVTimeKernel.C)This kernel implements the term:
The implementation is identical to the implementation of FEM kernels except that the FVM does not require multiplication by the test function.