ADKernels System

An "ADKernel" 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 ADKernel class.

In an ADKernel subclass the computeQpResidual() function must be overridden. This is where you implement your PDE weak form terms. By using an ADKernel as opposed to an ordinary Kernel, you are automatically computing the Jacobian, so you do not have to worry about hand-coding a Jacobian, which is prone to error.

Inside your ADKernel class, you have access to several member variables for computing the residual and Jacobian values in the above mentioned functions:

  • _i: indices for the current test shape functions.

  • _qp: current quadrature point index.

  • _u, _grad_u: value and gradient of the variable this ADKernel 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]).

  • _q_point: XYZ coordinates of the current quadrature point.

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

Custom ADKernel Creation

To create a custom kernel, you can follow the pattern of the ADDiffusion kernel implemented and included in the MOOSE framework.

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):

(2)

and then integrating by parts which gives the weak form:

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

The diffusion kernel header and implementation files are:










#ifndef ADDIFFUSION_H
#define ADDIFFUSION_H

#include "ADKernel.h"

template <ComputeStage compute_stage>
class ADDiffusion;

declareADValidParams(ADDiffusion);

template <ComputeStage compute_stage>
class ADDiffusion : public ADKernel<compute_stage>
{
public:
  ADDiffusion(const InputParameters & parameters);

protected:
  virtual ADResidual computeQpResidual() override;

  usingKernelMembers;
};

#endif /* ADDIFFUSION_H */
(framework/include/kernels/ADDiffusion.h)









#include "ADDiffusion.h"

registerADMooseObject("MooseApp", ADDiffusion);

defineADValidParams(
    ADDiffusion,
    ADKernel,
    params.addClassDescription("Same as `Diffusion` in terms of physics/residual, but the Jacobian "
                               "is computed using forward automatic differentiation"););

template <ComputeStage compute_stage>
ADDiffusion<compute_stage>::ADDiffusion(const InputParameters & parameters)
  : ADKernel<compute_stage>(parameters)
{
}

template <ComputeStage compute_stage>
ADResidual
ADDiffusion<compute_stage>::computeQpResidual()
{
  return _grad_u[_qp] * _grad_test[_i][_qp];
}
(framework/src/kernels/ADDiffusion.C)