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.

The Kernel system supports the use of AD for residual calculations, as such there are two options for creating Kernel objects: Kernel and ADKernel. To further understand automatic differentiation, please refer to the Automatic Differentiation page for more information.

In a Kernel subclass the computeQpResidual() function must be overridden. This is where you implement your PDE weak form terms. For non-AD objects 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.

Optimized Kernel Objects

Depending on the residual calculation being performed it is sometimes possible to optimize the calculation of the residual by precomputing values during the finite element assembly of the residual vector. The following table details the various Kernel base classes that can be used for as base classes to improve performance.

BaseOverrideUse
Kernel
ADKernel
computeQpResidualUse when the term in the PDE is multiplied by both the test function and the gradient of the test function (_test and _grad_test must be applied)
KernelValue
ADKernelValue
precomputeQpResidualUse when the term computed in the PDE is only multiplied by the test function (do not use _test in the override, it is applied automatically)
KernelGrad
ADKernelGrad
precomputeQpResidualUse when the term computed in the PDE is only multiplied by the gradient of the test function (do not use _grad_test in the override, it is applied automatically)

Custom Kernel Creation

To create a custom kernel, you can follow the pattern of the Diffusion or ADDiffusion objects implemented and included in the MOOSE framework. Additionally, Example 2 in MOOSE provides a step-by-step overview of creating your own custom kernel. The following describes that calculation of the diffusion term of a PDE.

The strong-form of the diffusion equation 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 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):

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)

As mentioned, the computeQpResidual method must be overridden for both flavors of kernels non-AD and AD. The computeQpResidual method for the non-AD version, Diffusion, is provided in Listing 1.

Listing 1: The C++ weak-form residual statement of Eq. (2) as implemented in the Diffusion kernel.

Real
Diffusion::computeQpResidual()
{
  return _grad_u[_qp] * _grad_test[_i][_qp];
}
(framework/src/kernels/Diffusion.C)

This object also overrides the computeQpJacobian method to define Jacobian term of (test/tests/test_harness/duplicate_outputs_analyzejacobian) as shown in Listing 2.

Listing 2: The C++ weak-form Jacobian statement of (test/tests/test_harness/duplicate_outputs_analyzejacobian) as implemented in the Diffusion kernel.

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

The AD version of this object, ADDiffusion, relies on an optimized kernel object (see Optimized Kernel Objects), as such it overrides precomputeQpResidual as follows.

Listing 3: The C++ pre-computed portions of the weak-form residual statement of Eq. (2) as implemented in the ADDiffusion kernel.

ADDiffusion::precomputeQpResidual()
{
  return _grad_u[_qp];
}
(framework/src/kernels/ADDiffusion.C)

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:

where is the finite element solution, and

(4)

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:

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

The derivative of equation Eq. (4) with respect to is then:

So that the Jacobian term for equation Eq. (4) is

where is what we call du_dot_du in MOOSE.

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


return _test[_i][_qp] * _u_dot[_qp];

And the corresponding computeQpJacobian() is:


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

Coupling with Scalar Variables

If the weak form has contributions from scalar variables, then this contribution can be treated similarly as coupling from other spatial variables. See the Coupleable interface for how to obtain the variable values. Residual contributions are simply added to the computeQpResidual() function. Jacobian terms from the test spatial variable and incremental scalar variable are added by overriding the function computeQpOffDiagJacobianScalar().

Contributions to the scalar variable weak equation (test scalar variable terms) are not natively treated by the Kernel class. Inclusion of these residual and Jacobian contributions are discussed within ScalarKernels and specifically KernelScalarBase.

Further Kernel Documentation

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

Available Objects

  • Moose App
  • ADBodyForceDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • ADCoefReactionImplements the residual term (p*u, test)
  • ADConservativeAdvectionConservative form of which in its weak form is given by: .
  • ADCoupledForceImplements a source term proportional to the value of a coupled variable. Weak form: .
  • ADCoupledTimeDerivativeTime derivative Kernel that acts on a coupled variable. Weak form: .
  • ADDiffusionSame as Diffusion in terms of physics/residual, but the Jacobian is computed using forward automatic differentiation
  • ADMatBodyForceKernel that defines a body force modified by a material property
  • ADMatCoupledForceKernel representing the contribution of the PDE term , where is a material property coefficient, is a coupled scalar field variable, and Jacobian derivatives are calculated using automatic differentiation.
  • ADMatDiffusionDiffusion equation kernel that takes an isotropic diffusivity from a material property
  • ADMatReactionKernel representing the contribution of the PDE term , where is a reaction rate material property, is a scalar variable (nonlinear or coupled), and whose Jacobian contribution is calculated using automatic differentiation.
  • ADMaterialPropertyValueResidual term (u - prop) to set variable u equal to a given material property prop
  • ADReactionImplements a simple consuming reaction term with weak form .
  • ADScalarLMKernelThis class is used to enforce integral of phi = V_0 with a Lagrange multiplier approach.
  • ADTimeDerivativeThe time derivative operator with the weak form of .
  • ADVectorDiffusionThe Laplacian operator (), with the weak form of . The Jacobian is computed using automatic differentiation
  • ADVectorTimeDerivativeThe time derivative operator with the weak form of .
  • AnisotropicDiffusionAnisotropic diffusion kernel with weak form given by .
  • ArrayBodyForceApplies body forces specified with functions to an array variable.
  • ArrayCoupledTimeDerivativeTime derivative Array Kernel that acts on a coupled variable. Weak form: . The coupled variable and the variable must have the same dimensionality
  • ArrayDiffusionThe array Laplacian operator (), with the weak form of .
  • ArrayReactionThe array reaction operator with the weak form of .
  • ArrayTimeDerivativeArray time derivative operator with the weak form of .
  • BodyForceDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • CoefReactionImplements the residual term (p*u, test)
  • 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 .
  • DivFieldThe divergence operator optionally scaled by a constant scalar coefficient. Weak form: .
  • FunctionDiffusionDiffusion with a function coefficient.
  • GradFieldThe gradient operator optionally scaled by a constant scalar coefficient. Weak form: .
  • 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 .
  • MassMatrixComputes a finite element mass matrix meant for use in preconditioning schemes which require one
  • MatBodyForceKernel that defines a body force modified by a material property
  • MatCoupledForceImplements a forcing term RHS of the form PDE = RHS, where RHS = Sum_j c_j * m_j * v_j. c_j, m_j, and v_j are provided as real coefficients, material properties, and coupled variables, respectively.
  • MatDiffusionDiffusion equation Kernel that takes an isotropic Diffusivity from a material property
  • MatReactionKernel to add -L*v, where L=reaction rate, v=variable
  • 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.
  • MaterialPropertyValueResidual term (u - prop) to set variable u equal to a given material property prop
  • NullKernelKernel that sets a zero residual.
  • ReactionImplements a simple consuming reaction term with weak form .
  • ScalarLMKernelThis class is used to enforce integral of phi = V_0 with a Lagrange multiplier approach.
  • ScalarLagrangeMultiplierThis class is used to enforce integral of phi = V_0 with a Lagrange multiplier approach.
  • 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 .
  • VectorCoupledTimeDerivativeTime derivative Kernel that acts on a coupled vector variable. Weak form: .
  • VectorDiffusionThe Laplacian operator (), with the weak form of .
  • VectorFunctionReactionKernel representing the contribution of the PDE term , where is a function coefficient and is a vector variable.
  • VectorTimeDerivativeThe time derivative operator with the weak form of .
  • Navier Stokes App
  • DistributedForceImplements a force term in the Navier Stokes momentum equation.
  • DistributedPowerImplements the power term of a specified force in the Navier Stokes energy equation.
  • INSADBoussinesqBodyForceComputes a body force for natural convection buoyancy.
  • INSADEnergyAdvectionThis class computes the residual and Jacobian contributions for temperature advection for a divergence free velocity field.
  • INSADEnergyAmbientConvectionComputes a heat source/sink due to convection from ambient surroundings.
  • INSADEnergyMeshAdvectionThis class computes the residual and Jacobian contributions for temperature advection from mesh velocity in an ALE simulation.
  • INSADEnergySUPGAdds the supg stabilization to the INS temperature/energy equation
  • INSADEnergySourceComputes an arbitrary volumetric heat source (or sink).
  • INSADGravityForceComputes a body force due to gravity.
  • INSADHeatConductionTimeDerivativeAD Time derivative term of the heat equation for quasi-constant specific heat and the density .
  • INSADMassThis class computes the mass equation residual and Jacobian contributions (the latter using automatic differentiation) for the incompressible Navier-Stokes equations.
  • INSADMassPSPGThis class adds PSPG stabilization to the mass equation, enabling use of equal order shape functions for pressure and velocity variables
  • INSADMomentumAdvectionAdds the advective term to the INS momentum equation
  • INSADMomentumCoupledForceComputes a body force due to a coupled vector variable or a vector function
  • INSADMomentumGradDivAdds grad-div stabilization to the INS momentum equation
  • INSADMomentumMeshAdvectionCorrects the convective derivative for situations in which the fluid mesh is dynamic.
  • INSADMomentumPressureAdds the pressure term to the INS momentum equation
  • INSADMomentumSUPGAdds the supg stabilization to the INS momentum equation
  • INSADMomentumTimeDerivativeThis class computes the time derivative for the incompressible Navier-Stokes momentum equation.
  • INSADMomentumViscousAdds the viscous term to the INS momentum equation
  • INSADSmagorinskyEddyViscosityComputes eddy viscosity term using Smagorinky's LES model
  • INSChorinCorrectorThis class computes the 'Chorin' Corrector equation in fully-discrete (both time and space) form.
  • INSChorinPredictorThis class computes the 'Chorin' Predictor equation in fully-discrete (both time and space) form.
  • INSChorinPressurePoissonThis class computes the pressure Poisson solve which is part of the 'split' scheme used for solving the incompressible Navier-Stokes equations.
  • INSCompressibilityPenaltyThe penalty term may be used when Dirichlet boundary condition is applied to the entire boundary.
  • INSFEFluidEnergyKernelAdds advection, diffusion, and heat source terms to energy equation, potentially with stabilization
  • INSFEFluidMassKernelAdds advective term of mass conservation equation along with pressure-stabilized Petrov-Galerkin terms
  • INSFEFluidMomentumKernelAdds advection, viscous, pressure, friction, and gravity terms to the Navier-Stokes momentum equation, potentially with stabilization
  • INSMassThis class computes the mass equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation.
  • INSMassRZThis class computes the mass equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation in RZ coordinates.
  • INSMomentumLaplaceFormThis class computes momentum equation residual and Jacobian viscous contributions for the 'Laplacian' form of the governing equations.
  • INSMomentumLaplaceFormRZThis class computes additional momentum equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation in RZ (axisymmetric cylindrical) coordinates, using the 'Laplace' form of the governing equations.
  • INSMomentumTimeDerivativeThis class computes the time derivative for the incompressible Navier-Stokes momentum equation.
  • INSMomentumTractionFormThis class computes momentum equation residual and Jacobian viscous contributions for the 'traction' form of the governing equations.
  • INSMomentumTractionFormRZThis class computes additional momentum equation residual and Jacobian contributions for the incompressible Navier-Stokes momentum equation in RZ (axisymmetric cylindrical) coordinates.
  • INSPressurePoissonThis class computes the pressure Poisson solve which is part of the 'split' scheme used for solving the incompressible Navier-Stokes equations.
  • INSProjectionThis class computes the 'projection' part of the 'split' method for solving incompressible Navier-Stokes.
  • INSSplitMomentumThis class computes the 'split' momentum equation residual.
  • INSTemperatureThis class computes the residual and Jacobian contributions for the incompressible Navier-Stokes temperature (energy) equation.
  • INSTemperatureTimeDerivativeThis class computes the time derivative for the incompressible Navier-Stokes momentum equation.
  • MDFluidEnergyKernelAdds advection, diffusion, and heat source terms to energy equation, potentially with stabilization
  • MDFluidMassKernelAdds advective term of mass conservation equation along with pressure-stabilized Petrov-Galerkin terms
  • MDFluidMomentumKernelAdds advection, viscous, pressure, friction, and gravity terms to the Navier-Stokes momentum equation, potentially with stabilization
  • MassConvectiveFluxImplements the advection term for the Navier Stokes mass equation.
  • MomentumConvectiveFluxImplements the advective term of the Navier Stokes momentum equation.
  • NSEnergyInviscidFluxThis class computes the inviscid part of the energy flux.
  • NSEnergyThermalFluxThis class is responsible for computing residuals and Jacobian terms for the k * grad(T) * grad(phi) term in the Navier-Stokes energy equation.
  • NSEnergyViscousFluxViscous flux terms in energy equation.
  • NSGravityForceThis class computes the gravity force contribution.
  • NSGravityPowerThis class computes the momentum contributed by gravity.
  • NSMassInviscidFluxThis class computes the inviscid flux in the mass equation.
  • NSMomentumInviscidFluxThe inviscid flux (convective + pressure terms) for the momentum conservation equations.
  • NSMomentumInviscidFluxWithGradPThis class computes the inviscid flux with pressure gradient in the momentum equation.
  • NSMomentumViscousFluxDerived instance of the NSViscousFluxBase class for the momentum equations.
  • NSSUPGEnergyCompute residual and Jacobian terms form the SUPG terms in the energy equation.
  • NSSUPGMassCompute residual and Jacobian terms form the SUPG terms in the mass equation.
  • NSSUPGMomentumCompute residual and Jacobian terms form the SUPG terms in the momentum equation.
  • NSTemperatureL2This class was originally used to solve for the temperature using an L2-projection.
  • PINSFEFluidPressureTimeDerivativeAdds the transient term of the porous-media mass conservation equation
  • PINSFEFluidTemperatureTimeDerivativeThe time derivative operator with the weak form of .
  • PINSFEFluidVelocityTimeDerivativeThe time derivative operator with the weak form of .
  • PMFluidPressureTimeDerivativeAdds the transient term of the porous-media mass conservation equation
  • PMFluidTemperatureTimeDerivativeThe time derivative operator with the weak form of .
  • PMFluidVelocityTimeDerivativeThe time derivative operator with the weak form of .
  • PressureGradientImplements the pressure gradient term for one of the Navier Stokes momentum equations.
  • TotalEnergyConvectiveFluxImplements the advection term for the Navier Stokes energy equation.
  • VectorMassMatrixComputes a finite element mass matrix meant for use in preconditioning schemes which require one
  • Misc App
  • ADThermoDiffusionCalculates diffusion due to temperature gradient and Soret Coefficient
  • CoefDiffusionKernel for diffusion with diffusivity = coef + function
  • ThermoDiffusionKernel for thermo-diffusion (Soret effect, thermophoresis, etc.)
  • XFEMApp
  • CrackTipEnrichmentStressDivergenceTensorsEnrich stress divergence kernel for small-strain simulations
  • Peridynamics App
  • ForceStabilizedSmallStrainMechanicsNOSPDClass for calculating the residual and Jacobian for the force-stabilized peridynamic correspondence model under small strain assumptions
  • GeneralizedPlaneStrainOffDiagNOSPDClass for calculating the off-diagonal Jacobian of the coupling between displacements (or temperature) with scalar out-of-plane strain for the generalized plane strain using the H1NOSPD formulation
  • GeneralizedPlaneStrainOffDiagOSPDClass for calculating the off-diagonal Jacobian corresponding to coupling between displacements (or temperature) and the scalar out-of-plane strain for the generalized plane strain using the OSPD formulation
  • HeatConductionBPDClass for calculating the residual and Jacobian for the bond-based peridynamic heat conduction formulation
  • HeatSourceBPDClass for calculating the residual from heat source for the bond-based peridynamic heat conduction formulation
  • HorizonStabilizedFormIFiniteStrainMechanicsNOSPDClass for calculating the residual and the Jacobian for Form I of the horizon-stabilized peridynamic correspondence model under finite strain assumptions
  • HorizonStabilizedFormIIFiniteStrainMechanicsNOSPDClass for calculating the residual and the Jacobian for Form II of the horizon-stabilized peridynamic correspondence model under finite strain assumptions
  • HorizonStabilizedFormIISmallStrainMechanicsNOSPDClass for calculating the residual and the Jacobian for Form II of the horizon-stabilized peridynamic correspondence model under small strain assumptions
  • HorizonStabilizedFormISmallStrainMechanicsNOSPDClass for calculating the residual and the Jacobian for Form I of the horizon-stabilizedperidynamic correspondence model under small strain assumptions
  • MechanicsBPDClass for calculating the residual and Jacobian for the bond-based peridynamic mechanics formulation
  • MechanicsOSPDClass for calculating the residual and Jacobian for the ordinary state-based peridynamic mechanics formulation
  • WeakPlaneStressNOSPDClass for calculating the residual and the Jacobian for the peridynamic plane stress model using weak formulation based on peridynamic correspondence models
  • Thermal Hydraulics App
  • ADHeatConductionRZAdds a heat conduction term in XY coordinates interpreted as cylindrical coordinates
  • ADHeatConductionTimeDerivativeRZAdds a time derivative term for the energy equation in XY coordinates interpreted as cylindrical coordinates
  • ADHeatStructureHeatSourceAdds a heat source term for the energy equation
  • ADHeatStructureHeatSourceRZAdds a heat source term in XY coordinates interpreted as cylindrical coordinates
  • ADOneD3EqnEnergyGravityComputes the gravity term for the energy equation in 1-phase flow
  • ADOneD3EqnEnergyHeatFluxComputes a heat flux term for the energy equation in a flow channel
  • ADOneD3EqnEnergyHeatFluxFromHeatStructure3DComputes a heat flux term from a 3D heat structure in the energy equation for 1-phase flow
  • ADOneD3EqnMomentumAreaGradientComputes the area gradient term in the momentum equation for single phase flow.
  • ADOneD3EqnMomentumFormLossComputes a volumetric form loss for the momentum equation for 1-phase flow
  • ADOneD3EqnMomentumFrictionComputes wall friction term for single phase flow.
  • ADOneD3EqnMomentumGravityComputes gravity term for the momentum equation for 1-phase flow
  • ADOneDEnergyWallHeatFluxComputes a heat flux term for the energy equation
  • ADOneDEnergyWallHeatingComputes a convective heat flux term for the energy equation for 1-phase flow
  • CoupledForceRZAdds a coupled force term in XY coordinates interpreted as cylindrical coordinates
  • OneD3EqnEnergyFluxComputes an energy flux for single phase flow
  • OneD3EqnEnergyGravityComputes a gravity term for the energy equation in 1-phase flow
  • OneD3EqnEnergyHeatSourceComputes a volumetric heat source for 1-phase flow channel
  • OneD3EqnMomentumAreaGradientComputes the area gradient term in the momentum equation for single phase flow.
  • OneD3EqnMomentumFluxComputes a momentum flux term for 1-phase flow
  • OneD3EqnMomentumFormLossComputes a form loss term for the momentum equation for 1-phase flow
  • OneD3EqnMomentumFrictionComputes wall friction term for single phase flow.
  • OneD3EqnMomentumGravityComputes gravity term for the momentum equation for 1-phase flow
  • OneDEnergyWallHeatFluxAdds a heat flux along the local heated perimeter
  • OneDEnergyWallHeatingAdds a convective heat flux term from a wall temperature
  • Solid Mechanics App
  • ADDynamicStressDivergenceTensorsResidual due to stress related Rayleigh damping and HHT time integration terms
  • ADGravityApply gravity. Value is in units of acceleration.
  • ADInertialForceCalculates the residual for the inertial force () and the contribution of mass dependent Rayleigh damping and HHT time integration scheme ($\eta \cdot M \cdot ((1+\alpha)velq2-\alpha \cdot vel-old) $)
  • ADInertialForceShellCalculates the residual for the inertial force/moment and the contribution of mass dependent Rayleigh damping and HHT time integration scheme.
  • ADStressDivergenceRSphericalTensorsCalculate stress divergence for a spherically symmetric 1D problem in polar coordinates.
  • ADStressDivergenceRZTensorsCalculate stress divergence for an axisymmetric problem in cylindrical coordinates.
  • ADStressDivergenceShellQuasi-static stress divergence kernel for Shell element
  • ADStressDivergenceTensorsStress divergence kernel with automatic differentiation for the Cartesian coordinate system
  • ADSymmetricStressDivergenceTensorsStress divergence kernel with automatic differentiation for the Cartesian coordinate system
  • ADWeakPlaneStressPlane stress kernel to provide out-of-plane strain contribution.
  • AsymptoticExpansionHomogenizationKernelKernel for asymptotic expansion homogenization for elasticity
  • CosseratStressDivergenceTensorsStress divergence kernel for the Cartesian coordinate system
  • DynamicStressDivergenceTensorsResidual due to stress related Rayleigh damping and HHT time integration terms
  • GeneralizedPlaneStrainOffDiagGeneralized Plane Strain kernel to provide contribution of the out-of-plane strain to other kernels
  • GravityApply gravity. Value is in units of acceleration.
  • HomogenizedTotalLagrangianStressDivergenceTotal Lagrangian stress equilibrium kernel with homogenization constraint Jacobian terms
  • InertialForceCalculates the residual for the inertial force () and the contribution of mass dependent Rayleigh damping and HHT time integration scheme ($\eta \cdot M \cdot ((1+\alpha)velq2-\alpha \cdot vel-old) $)
  • InertialForceBeamCalculates the residual for the inertial force/moment and the contribution of mass dependent Rayleigh damping and HHT time integration scheme.
  • InertialTorqueKernel for inertial torque: density * displacement x acceleration
  • MaterialVectorBodyForceApply a body force vector to the coupled displacement component.
  • MomentBalancing
  • OutOfPlanePressureApply pressure in the out-of-plane direction in 2D plane stress or generalized plane strain models
  • PhaseFieldFractureMechanicsOffDiagStress divergence kernel for phase-field fracture: Computes off diagonal damage dependent Jacobian components. To be used with StressDivergenceTensors or DynamicStressDivergenceTensors.
  • PlasticHeatEnergyPlastic heat energy density = coeff * stress * plastic_strain_rate
  • PoroMechanicsCouplingAdds , where the subscript is the component.
  • StressDivergenceBeamQuasi-static and dynamic stress divergence kernel for Beam element
  • StressDivergenceRSphericalTensorsCalculate stress divergence for a spherically symmetric 1D problem in polar coordinates.
  • StressDivergenceRZTensorsCalculate stress divergence for an axisymmetric problem in cylindrical coordinates.
  • StressDivergenceTensorsStress divergence kernel for the Cartesian coordinate system
  • StressDivergenceTensorsTrussKernel for truss element
  • TotalLagrangianStressDivergenceEnforce equilibrium with a total Lagrangian formulation in Cartesian coordinates.
  • TotalLagrangianStressDivergenceAxisymmetricCylindricalEnforce equilibrium with a total Lagrangian formulation in axisymmetric cylindrical coordinates.
  • TotalLagrangianStressDivergenceCentrosymmetricSphericalEnforce equilibrium with a total Lagrangian formulation in centrosymmetric spherical coordinates.
  • TotalLagrangianWeakPlaneStressPlane stress kernel to provide out-of-plane strain contribution.
  • UpdatedLagrangianStressDivergenceEnforce equilibrium with an updated Lagrangian formulation in Cartesian coordinates.
  • WeakPlaneStressPlane stress kernel to provide out-of-plane strain contribution.
  • Electromagnetics App
  • CurlCurlFieldWeak form term corresponding to .
  • VectorCurrentSourceKernel to calculate the current source term in the Helmholtz wave equation.
  • VectorSecondTimeDerivativeThe second time derivative operator for vector variables.
  • Porous Flow App
  • FluxLimitedTVDAdvectionConservative form of (advection), using the Flux Limited TVD scheme invented by Kuzmin and Turek
  • PorousFlowAdvectiveFluxFully-upwinded advective flux of the component given by fluid_component
  • PorousFlowBasicAdvectionAdvective flux of a Variable using the Darcy velocity of the fluid phase
  • PorousFlowDesorpedMassTimeDerivativeDesorped component mass derivative wrt time.
  • PorousFlowDesorpedMassVolumetricExpansionDesorped_mass * rate_of_solid_volumetric_expansion
  • PorousFlowDispersiveFluxDispersive and diffusive flux of the component given by fluid_component in all phases
  • PorousFlowEffectiveStressCouplingImplements the weak form of the expression biot_coefficient * grad(effective fluid pressure)
  • PorousFlowEnergyTimeDerivativeDerivative of heat-energy-density wrt time
  • PorousFlowExponentialDecayResidual = rate * (variable - reference). Useful for modelling exponential decay of a variable
  • PorousFlowFluxLimitedTVDAdvectionAdvective flux of fluid species or heat using the Flux Limited TVD scheme invented by Kuzmin and Turek
  • PorousFlowFullySaturatedAdvectiveFluxFully-upwinded advective flux of the fluid component given by fluid_component, in a single-phase fluid
  • PorousFlowFullySaturatedDarcyBaseDarcy flux suitable for models involving a fully-saturated, single phase, single component fluid. No upwinding is used
  • PorousFlowFullySaturatedDarcyFlowDarcy flux suitable for models involving a fully-saturated single phase, multi-component fluid. No upwinding is used
  • PorousFlowFullySaturatedHeatAdvectionHeat flux that arises from the advection of a fully-saturated single phase fluid. No upwinding is used
  • PorousFlowFullySaturatedMassTimeDerivativeFully-saturated version of the single-component, single-phase fluid mass derivative wrt time
  • PorousFlowFullySaturatedUpwindHeatAdvectionHeat advection by a fluid. The fluid is assumed to have a single phase, and the advection is fully upwinded
  • PorousFlowHeatAdvectionFully-upwinded heat flux, advected by the fluid
  • PorousFlowHeatConductionHeat conduction in the Porous Flow module
  • PorousFlowHeatMassTransferCalculate heat or mass transfer from a coupled variable v to the variable u. No mass lumping is performed here.
  • PorousFlowHeatVolumetricExpansionEnergy-density*rate_of_solid_volumetric_expansion. The energy-density is lumped to the nodes
  • PorousFlowMassRadioactiveDecayRadioactive decay of a fluid component
  • PorousFlowMassTimeDerivativeDerivative of fluid-component mass with respect to time. Mass lumping to the nodes is used.
  • PorousFlowMassVolumetricExpansionComponent_mass*rate_of_solid_volumetric_expansion. This Kernel lumps the component mass to the nodes.
  • PorousFlowPlasticHeatEnergyPlastic heat energy density source = (1 - porosity) * coeff * stress * plastic_strain_rate
  • PorousFlowPreDisPrecipitation-dissolution of chemical species
  • Phase Field App
  • ACBarrierFunctionAllen-Cahn kernel used when 'mu' is a function of variables
  • ACGBPolyGrain-Boundary model concentration dependent residual
  • ACGrGrElasticDrivingForceAdds elastic energy contribution to the Allen-Cahn equation
  • ACGrGrMultiMulti-phase poly-crystalline Allen-Cahn Kernel
  • ACGrGrPolyGrain-Boundary model poly-crystalline interface Allen-Cahn Kernel
  • ACGrGrPolyLinearizedInterfaceGrain growth model Allen-Cahn Kernel with linearized interface variable transformation
  • ACInterfaceGradient energy Allen-Cahn Kernel
  • ACInterface2DMultiPhase1Gradient energy Allen-Cahn Kernel where the derivative of interface parameter kappa wrt the gradient of order parameter is considered.
  • ACInterface2DMultiPhase2Gradient energy Allen-Cahn Kernel where the interface parameter kappa is considered.
  • ACInterfaceChangedVariableGradient energy Allen-Cahn Kernel using a change of variable
  • ACInterfaceCleavageFractureGradient energy Allen-Cahn Kernel where crack propagation along weakcleavage plane is preferred
  • ACInterfaceKobayashi1Anisotropic gradient energy Allen-Cahn Kernel Part 1
  • ACInterfaceKobayashi2Anisotropic Gradient energy Allen-Cahn Kernel Part 2
  • ACInterfaceStressInterface stress driving force Allen-Cahn Kernel
  • ACKappaFunctionGradient energy term for when kappa as a function of the variable
  • ACMultiInterfaceGradient energy Allen-Cahn Kernel with cross terms
  • ACSEDGPolyStored Energy contribution to grain growth
  • ACSwitchingKernel for Allen-Cahn equation that adds derivatives of switching functions and energies
  • ADACBarrierFunctionAllen-Cahn kernel used when 'mu' is a function of variables
  • ADACGrGrMultiMulti-phase poly-crystalline Allen-Cahn Kernel
  • ADACInterfaceGradient energy Allen-Cahn Kernel
  • ADACInterfaceKobayashi1Anisotropic gradient energy Allen-Cahn Kernel Part 1
  • ADACInterfaceKobayashi2Anisotropic Gradient energy Allen-Cahn Kernel Part 2
  • ADACKappaFunctionGradient energy term for when kappa as a function of the variable
  • ADACSwitchingKernel for Allen-Cahn equation that adds derivatives of switching functions and energies
  • ADAllenCahnAllen-Cahn Kernel that uses a DerivativeMaterial Free Energy
  • ADCHSoretMobilityAdds contribution due to thermo-migration to the Cahn-Hilliard equation using a concentration 'u', temperature 'T', and thermal mobility 'mobility' (in units of length squared per time).
  • ADCHSplitChemicalPotentialChemical potential kernel in Split Cahn-Hilliard that solves chemical potential in a weak form
  • ADCHSplitConcentrationConcentration kernel in Split Cahn-Hilliard that solves chemical potential in a weak form
  • ADCoefCoupledTimeDerivativeScaled time derivative Kernel that acts on a coupled variable
  • ADCoupledSwitchingTimeDerivativeCoupled time derivative Kernel that multiplies the time derivative by
  • ADGrainGrowthGrain-Boundary model poly-crystalline interface Allen-Cahn Kernel
  • ADMatAnisoDiffusionDiffusion equation kernel that takes an anisotropic diffusivity from a material property
  • ADSplitCHParsedSplit formulation Cahn-Hilliard Kernel that uses a DerivativeMaterial Free Energy
  • ADSplitCHWResSplit formulation Cahn-Hilliard Kernel for the chemical potential variable with a scalar (isotropic) mobility
  • ADSplitCHWResAnisoSplit formulation Cahn-Hilliard Kernel for the chemical potential variable with a scalar (isotropic) mobility
  • ADSusceptibilityTimeDerivativeA modified time derivative Kernel that multiplies the time derivative of a variable by a generalized susceptibility
  • AllenCahnAllen-Cahn Kernel that uses a DerivativeMaterial Free Energy
  • AllenCahnElasticEnergyOffDiagThis kernel calculates off-diagonal Jacobian of elastic energy in AllenCahn with respect to displacements
  • AntitrappingCurrentKernel that provides antitrapping current at the interface for alloy solidification
  • CHBulkPFCTradCahn-Hilliard kernel for a polynomial phase field crystal free energy.
  • CHInterfaceGradient energy Cahn-Hilliard Kernel with a scalar (isotropic) mobility
  • CHInterfaceAnisoGradient energy Cahn-Hilliard Kernel with a tensor (anisotropic) mobility
  • CHMathSimple demonstration Cahn-Hilliard Kernel using an algebraic double-well potential
  • CHPFCRFFCahn-Hilliard residual for the RFF form of the phase field crystal model
  • CHSplitChemicalPotentialChemical potential kernel in Split Cahn-Hilliard that solves chemical potential in a weak form
  • CHSplitConcentrationConcentration kernel in Split Cahn-Hilliard that solves chemical potential in a weak form
  • CHSplitFluxComputes flux as nodal variable
  • CahnHilliardCahn-Hilliard Kernel that uses a DerivativeMaterial Free Energy and a scalar (isotropic) mobility
  • CahnHilliardAnisoCahn-Hilliard Kernel that uses a DerivativeMaterial Free Energy and a tensor (anisotropic) mobility
  • ChangedVariableTimeDerivativeA modified time derivative Kernel that multiplies the time derivative bythe derivative of the nonlinear preconditioning function
  • CoefCoupledTimeDerivativeScaled time derivative Kernel that acts on a coupled variable
  • ConservedLangevinNoiseSource term for noise from a ConservedNoise userobject
  • CoupledAllenCahnCoupled Allen-Cahn Kernel that uses a DerivativeMaterial Free Energy
  • CoupledMaterialDerivativeKernel that implements the first derivative of a function material property with respect to a coupled variable.
  • CoupledSusceptibilityTimeDerivativeA modified coupled time derivative Kernel that multiplies the time derivative of a coupled variable by a generalized susceptibility
  • CoupledSwitchingTimeDerivativeCoupled time derivative Kernel that multiplies the time derivative by
  • DiscreteNucleationForceTerm for inserting grain nuclei or phases in non-conserved order parameter fields
  • GradientComponentSet the kernel variable to a specified component of the gradient of a coupled variable.
  • HHPFCRFFReaction type kernel for the RFF phase fit crystal model
  • KKSACBulkCKKS model kernel (part 2 of 2) for the Bulk Allen-Cahn. This includes all terms dependent on chemical potential.
  • KKSACBulkFKKS model kernel (part 1 of 2) for the Bulk Allen-Cahn. This includes all terms NOT dependent on chemical potential.
  • KKSCHBulkKKS model kernel for the Bulk Cahn-Hilliard term. This operates on the concentration 'c' as the non-linear variable
  • KKSMultiACBulkCMulti-phase KKS model kernel (part 2 of 2) for the Bulk Allen-Cahn. This includes all terms dependent on chemical potential.
  • KKSMultiACBulkFKKS model kernel (part 1 of 2) for the Bulk Allen-Cahn. This includes all terms NOT dependent on chemical potential.
  • KKSMultiPhaseConcentrationKKS multi-phase model kernel to enforce . The non-linear variable of this kernel is , the final phase concentration in the list.
  • KKSPhaseChemicalPotentialKKS model kernel to enforce the pointwise equality of phase chemical potentials . The non-linear variable of this kernel is .
  • KKSPhaseConcentrationKKS model kernel to enforce the decomposition of concentration into phase concentration . The non-linear variable of this kernel is .
  • KKSSplitCHCResKKS model kernel for the split Bulk Cahn-Hilliard term. This kernel operates on the physical concentration 'c' as the non-linear variable
  • LangevinNoiseSource term for non-conserved Langevin noise
  • LaplacianSplitSplit with a variable that holds the Laplacian of a phase field variable.
  • MaskedBodyForceKernel that defines a body force modified by a material mask
  • MaskedExponentialKernel to add dilute solution term to Poisson's equation for electrochemical sintering
  • MatAnisoDiffusionDiffusion equation Kernel that takes an anisotropic Diffusivity from a material property
  • MatGradSquareCoupledGradient square of a coupled variable.
  • MultiGrainRigidBodyMotionAdds rigid body motion to grains
  • NestedKKSACBulkCKKS model kernel (part 2 of 2) for the Bulk Allen-Cahn. This includes all terms dependent on chemical potential.
  • NestedKKSACBulkFKKS model kernel (part 1 of 2) for the Bulk Allen-Cahn. This includes all terms NOT dependent on chemical potential.
  • NestedKKSSplitCHCResKKS model kernel for the split Bulk Cahn-Hilliard term. This kernel operates on the physical concentration 'c' as the non-linear variable.
  • SLKKSChemicalPotentialSLKKS model kernel to enforce the pointwise equality of sublattice chemical potentials in the same phase.
  • SLKKSMultiACBulkCMulti-phase SLKKS model kernel for the bulk Allen-Cahn. This includes all terms dependent on chemical potential.
  • SLKKSMultiPhaseConcentrationSLKKS multi-phase model kernel to enforce . The non-linear variable of this kernel is a phase's sublattice concentration
  • SLKKSPhaseConcentrationSublattice KKS model kernel to enforce the decomposition of concentration into phase and sublattice concentrations The non-linear variable of this kernel is a sublattice concentration of phase b.
  • SLKKSSumEnforce the sum of sublattice concentrations to a given phase concentration.
  • SimpleACInterfaceGradient energy for Allen-Cahn Kernel with constant Mobility and Interfacial parameter
  • SimpleCHInterfaceGradient energy for Cahn-Hilliard equation with constant Mobility and Interfacial parameter
  • SimpleCoupledACInterfaceGradient energy for Allen-Cahn Kernel with constant Mobility and Interfacial parameter for a coupled order parameter variable.
  • SimpleSplitCHWResGradient energy for split Cahn-Hilliard equation with constant Mobility for a coupled order parameter variable.
  • SingleGrainRigidBodyMotionAdds rigid mody motion to a single grain
  • SoretDiffusionAdd Soret effect to Split formulation Cahn-Hilliard Kernel
  • SplitCHMathSimple demonstration split formulation Cahn-Hilliard Kernel using an algebraic double-well potential
  • SplitCHParsedSplit formulation Cahn-Hilliard Kernel that uses a DerivativeMaterial Free Energy
  • SplitCHWResSplit formulation Cahn-Hilliard Kernel for the chemical potential variable with a scalar (isotropic) mobility
  • SplitCHWResAnisoSplit formulation Cahn-Hilliard Kernel for the chemical potential variable with a tensor (anisotropic) mobility
  • SusceptibilityTimeDerivativeA modified time derivative Kernel that multiplies the time derivative of a variable by a generalized susceptibility
  • SwitchingFunctionConstraintEtaLagrange multiplier kernel to constrain the sum of all switching functions in a multiphase system. This kernel acts on a non-conserved order parameter eta_i.
  • SwitchingFunctionConstraintLagrangeLagrange multiplier kernel to constrain the sum of all switching functions in a multiphase system. This kernel acts on the Lagrange multiplier variable.
  • SwitchingFunctionPenaltyPenalty kernel to constrain the sum of all switching functions in a multiphase system.
  • Heat Transfer App
  • ADHeatConductionSame as Diffusion in terms of physics/residual, but the Jacobian is computed using forward automatic differentiation
  • ADHeatConductionTimeDerivativeAD Time derivative term of the heat equation for quasi-constant specific heat and the density .
  • ADJouleHeatingSourceCalculates the heat source term corresponding to electrostatic Joule heating, with Jacobian contributions calculated using the automatic differentiation system.
  • ADMatHeatSourceForce term in thermal transport to represent a heat source
  • AnisoHeatConductionAnisotropic HeatConduction kernel with weak form given by .
  • AnisoHomogenizedHeatConductionKernel for asymptotic expansion homogenization for thermal conductivity when anisotropic thermal conductivities are used
  • HeatCapacityConductionTimeDerivativeTime derivative term of the heat equation with the heat capacity as an argument.
  • HeatConductionComputes residual/Jacobian contribution for term.
  • HeatConductionTimeDerivativeTime derivative term of the heat equation for quasi-constant specific heat and the density .
  • HeatSourceDemonstrates the multiple ways that scalar values can be introduced into kernels, e.g. (controllable) constants, functions, and postprocessors. Implements the weak form .
  • HomogenizedHeatConductionKernel for asymptotic expansion homogenization for thermal conductivity
  • JouleHeatingSourceCalculates the heat source term corresponding to electrostatic Joule heating.
  • SpecificHeatConductionTimeDerivativeTime derivative term of the heat equation with the specific heat and the density as arguments.
  • TrussHeatConductionComputes conduction term in heat equation for truss elements, taking cross-sectional area into account
  • TrussHeatConductionTimeDerivativeComputes time derivative term in heat equation for truss elements, taking cross-sectional area into account
  • Scalar Transport App
  • BodyForceLMImposes a body force onto a Lagrange multiplier constrained primal equation
  • CoupledForceLMAdds a coupled force term to a Lagrange multiplier constrained primal equation
  • LMDiffusionAdds a diffusion term to a Lagrange multiplier constrained primal equation
  • TimeDerivativeLMAdds a time derivative term to a Lagrange multiplier constrained primal equation
  • Chemical Reactions App
  • CoupledBEEquilibriumSubDerivative of equilibrium species concentration wrt time
  • CoupledBEKineticDerivative of kinetic species concentration wrt time
  • CoupledConvectionReactionSubConvection of equilibrium species
  • CoupledDiffusionReactionSubDiffusion of equilibrium species
  • DarcyFluxPressureDarcy flux: - cond * (Grad P - rho * g) where cond is the hydraulic conductivity, P is fluid pressure, rho is fluid density and g is gravity
  • DesorptionFromMatrixMass flow rate from the matrix to the porespace. Add this to TimeDerivative kernel to get complete DE for the fluid adsorbed in the matrix
  • DesorptionToPorespaceMass flow rate to the porespace from the matrix. Add this to the other kernels for the porepressure variable to form the complete DE
  • PrimaryConvectionConvection of primary species
  • PrimaryDiffusionDiffusion of primary species
  • PrimaryTimeDerivativeDerivative of primary species concentration wrt time
  • Geochemistry App
  • GeochemistryDispersionKernel describing grad(porosity * tensor_coeff * grad(concentration)), where porosity is an AuxVariable (or just a real number), tensor_coeff is the hydrodynamic dispersion tensor and concentration is the 'variable' for this Kernel
  • GeochemistryTimeDerivativeKernel describing porosity * d(concentration)/dt, where porosity is an AuxVariable (or just a real number) and concentration is the 'variable' for this Kernel. This Kernel should not be used if porosity is time-dependent. Mass lumping is employed for numerical stability
  • Fsi App
  • AcousticInertiaCalculates the residual for the inertial force which is the double time derivative of pressure.
  • ConvectedMeshCorrects the convective derivative for situations in which the fluid mesh is dynamic.
  • ConvectedMeshPSPGCorrects the convective derivative for situations in which the fluid mesh is dynamic.
  • Level Set App
  • LevelSetAdvectionImplements the level set advection equation: , where the weak form is .
  • LevelSetAdvectionSUPGSUPG stablization term for the advection portion of the level set equation.
  • LevelSetForcingFunctionSUPGThe SUPG stablization term for a forcing function.
  • LevelSetOlssonReinitializationThe re-initialization equation defined by Olsson et. al. (2007).
  • LevelSetTimeDerivativeSUPGSUPG stablization terms for the time derivative of the level set equation.