Kernels System

A "Kernel" is a piece of physics. It can represent one or more operators or terms in the weak form of a partial differential equation. With all terms on the left-hand-side, their sum is referred to as the "residual". The residual is evaluated at several integration quadrature points over the problem domain. To implement your own physics in MOOSE, you create your own kernel by subclassing the MOOSE Kernel class.

In a Kernel subclass the computeQpResidual() function must be overridden. This is where you implement your PDE weak form terms. The following member functions can optionally be overridden:

  • computeQpJacobian()

  • computeQpOffDiagJacobian()

These two functions provide extra information that can help the numerical solver(s) converge faster and better. Inside your Kernel class, you have access to several member variables for computing the residual and Jacobian values in the above mentioned functions:

  • _i, _j: indices for the current test and trial shape functions respectively.

  • _qp: current quadrature point index.

  • _u, _grad_u: value and gradient of the variable this Kernel operates on; indexed by _qp (i.e. _u[_qp]).

  • _test, _grad_test: value () and gradient () of the test functions at the q-points; indexed by _i and then _qp (i.e., _test[_i][_qp]).

  • _phi, _grad_phi: value () and gradient () of the trial functions at the q-points; indexed by _j and then _qp (i.e., _phi[_j][_qp]).

  • _q_point: XYZ coordinates of the current quadrature point.

  • _current_elem: pointer to the current element being operated on.

Custom Kernel Creation

To create a custom kernel, you can follow the pattern of the Diffusion kernel implemented and included in the MOOSE framework. Additionally, Example 2 provides a step-by-step overview of creating your own custom kernel.

The strong-form of the diffusion equations is defined on a 3-D domain as: find such that

(1) where is defined as the boundary on which the value of is fixed to a known constant , is defined as the boundary on which the flux across the boundary is fixed to a known constant , and $\hat{n} is the boundary outward normal.

The weak form is generated by multiplying by a test function () and integrating over the domain (using inner-product notation):

(-\nabla\cdot\nabla u_h, \psi_i) - (f, \psi_i) = 0\quad \forall\,\psi_i and then integrating by parts which gives the weak form:

(2) where is known as the trial function that defines the finite element discretization, , with being the basis functions.

The Jacobian, which is the derivative of Eq. 2 with respect to , is defined as:

(3)

The diffusion kernel header and implementation files are:

// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html

#ifndef DIFFUSION_H
#define DIFFUSION_H

#include "Kernel.h"

class Diffusion;

template <>
InputParameters validParams<Diffusion>();

/**
 * This kernel implements the Laplacian operator:
 * $\nabla u \cdot \nabla \phi_i$
 */
class Diffusion : public Kernel
{
public:
  Diffusion(const InputParameters & parameters);

protected:
  virtual Real computeQpResidual() override;

  virtual Real computeQpJacobian() override;
};

#endif /* DIFFUSION_H */
(../moose/framework/include/kernels/Diffusion.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html

#include "Diffusion.h"

registerMooseObject("MooseApp", Diffusion);

template <>
InputParameters
validParams<Diffusion>()
{
  InputParameters params = validParams<Kernel>();
  params.addClassDescription("The Laplacian operator ($-\\nabla \\cdot \\nabla u$), with the weak "
                             "form of $(\\nabla \\phi_i, \\nabla u_h)$.");
  return params;
}

Diffusion::Diffusion(const InputParameters & parameters) : Kernel(parameters) {}

Real
Diffusion::computeQpResidual()
{
  return _grad_u[_qp] * _grad_test[_i][_qp];
}

Real
Diffusion::computeQpJacobian()
{
  return _grad_phi[_j][_qp] * _grad_test[_i][_qp];
}
(../moose/framework/src/kernels/Diffusion.C)

Before a custom physics kernel is available for use, it must be registered in your application. This is done in e.g. src/base/YourApp.C for the YourApp application.


// src/base/YourApp.C contents

#include "YourKernel.h"

...

int YourApp::registerObjects(Factory & factory)
{
  ...
  registerKernel(YourKernel);
  ...
}

...

Time Derivative Kernels

You can create a time-derivative term/kernel by subclassing TimeKernel instead of Kernel. For example, the residual contribution for a time derivative term is:

(4)

where is the finite element solution, and

(5)

because you can interchange the order of differentiation and summation.

In the equation above, is the time derivative of the th finite element coefficient of . While the exact form of this derivative depends on the time stepping scheme, without much loss of generality, we can assume the following form for the time derivative:

(6)

for some constants , which depend on and the timestepping method.

The derivative of equation Eq. 5 with respect to is then:

(7)

So that the Jacobian term for equation Eq. 5 is

(8)

where is what we call du_dot_du in MOOSE.

Therefore the computeQpResidual() function for our time-derivative term kernel looks like:

text cpp return _test[_i][_qp] * _u_dot[_qp]; ``

And the corresponding computeQpJacobian() is:


return _test[_i][_qp] * _phi[_j][_qp] * _du_dot_du[_qp];

Further Kernel Documentation

Several specialized kernel types exist in MOOSE each with useful functionality. Details for each are in the sections below.

  • NumbatConvectionConvection of concentration with velocity given by Darcy's law
  • NumbatConvectionSFConvection of concentration with velocity given by Darcy's law using the streamfunction formulation
  • NumbatDarcyDarcy's law
  • NumbatDarcySFDarcy's law for the streamfunction formulation
  • NumbatDiffusionDiffusion kernel with porosity
  • NumbatDiffusionSFDiffusion kernel for the streamfunction formulation
  • NumbatTimeDerivativeTime derivative kernel
  • AnisotropicDiffusionAnisotropic diffusion kernel with weak form given by .
  • BodyForceDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • CoefTimeDerivativeThe time derivative operator with the weak form of .
  • ConservativeAdvectionConservative form of which in its weak form is given by: .
  • CoupledForceImplements a source term proportional to the value of a coupled variable. Weak form: .
  • CoupledTimeDerivativeTime derivative Kernel that acts on a coupled variable. Weak form: .
  • DiffusionThe Laplacian operator (), with the weak form of .
  • MassEigenKernelAn eigenkernel with weak form where is the eigenvalue.
  • MassLumpedTimeDerivativeLumped formulation of the time derivative . Its corresponding weak form is where denotes the time derivative of the solution coefficient associated with node .
  • MaterialDerivativeRankFourTestKernelClass used for testing derivatives of a rank four tensor material property.
  • MaterialDerivativeRankTwoTestKernelClass used for testing derivatives of a rank two tensor material property.
  • MaterialDerivativeTestKernelClass used for testing derivatives of a scalar material property.
  • NullKernelKernel that sets a zero residual.
  • ReactionImplements a simple consuming reaction term with weak form .
  • TimeDerivativeThe time derivative operator with the weak form of .
  • UserForcingFunctionDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • VectorBodyForceDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • VectorDiffusionThe Laplacian operator (), with the weak form of .