AuxKernels System

The AuxKernel system mimics the Kernels System but compute values that can be defined explicitly with a known function. There are two main use cases for the AuxKernel system: computing a quantity that varies with space and time for postprocessing purposes or for decoupling systems of equations. Examples for both of these use cases shall be discussed further in the following sections.

Creating a custom AuxKernel object is done by creating a new C++ object that inherits from AuxKernel, VectorAuxKernel or ArrayAuxKernel and overriding the computeValue method, which returns a scalar (Real), vector (RealVectorValue) or a Eigen vector (RealEigenVector) for the three types respectively. A forth type (AuxScalarKernel) also exists, but the syntax for these objects is different and detailed in the AuxScalarKernels System.

AuxKernel objects, like Kernel objects, must operate on a variable. Thus, there is a required parameter ("variable") that indicates the variable that the AuxKernel object is computing. These variables are defined in the AuxVariables block of the input file. AuxKernel objects derived from AuxKernel, VectorAuxKernel or ArrayAuxKernel operate on standard scalar, vector or array field variables respectively. For example the following input file snippet creates an auxiliary variable suitable for use with an VectorAuxKernel.

[AuxVariables]
  [vec]
    family = LAGRANGE_VEC
    order = FIRST
  []
[]
(test/tests/auxkernels/vector_function_aux/vector_function_aux.i)

Nodal vs Elemental AuxKernel Objects

There are two flavors of AuxKernel objects: nodal and elemental. The distinction is based on the type of variable that is being operated on by the object. If the variable family is LAGRANGE or LAGRANGE_VEC then the AuxKernel will behave as nodal. If the variable family is MONOMIAL then the AuxKernel will behave as elemental.

The difference is based on how the computeValue method of the object is called when the kernel is executed. In the nodal case the computeValue method will be executed on each node within the finite element mesh and the value returned from the method will directly assign the value of the shape function at that node.

In the elemental case the computeValue method will be executed on each quadrature point of an element of the finite element mesh. The values computed at the quadrature points are used to perform the correct finite element interpolation automatically and set the values for the degrees of freedom. Typically, in the elemental case the order of the monomial finite element is set to constant so there is a single DOF per element, but higher monomials are also supported.

As is evident by the functionality detailed, the distinction between the two arises from the nature of the finite element shape functions. For Lagrange shape functions the DOF values correspond with the nodes, while for elemental shape functions the DOF values are not associated with nodes.

The same AuxKernel object can be designed work both as elemental or nodal, for example the computeValue method for the FunctionAux object properly handles using the correct spatial location based on if the object is nodal or elemental with the isNodal method.

Block vs Boundary Restricted AuxKernel Objects

While auxiliary variables are always defined on mesh subdomains, MOOSE allows auxiliary kernels to be either block (mesh subdomain) or boundary restricted. When an auxiliary kernel is boundary restricted, it evaluates an auxiliary variable only on the designated boundaries. Because of this, the auxiliary variable will only have meaningful values on the boundaries even though it is defined on mesh subdomains. When an auxiliary kernel is block restricted, the variable that it evaluates must be defined on a subdomain covering the blocks where the auxiliary kernel is defined. When an auxiliary kernel is boundary restricted, the variable must be defined on a subdomain that all the sides on the boundaries are connected with. An elemental auxiliary variable defined on an element that has multiple boundary sides cannot be properly evaluated within a boundary restricted auxiliary kernel because elemental auxiliary variables can only store one value per element. Users can split the boundaries and define multiple elemental auxiliary variables for each split to avoid the situation of element connecting with multiple boundary sides.

Real
FunctionAux::computeValue()
{
  if (isNodal())
    return _func.value(_t, *_current_node);
  else
    return _func.value(_t, _q_point[_qp]);
}
(framework/src/auxkernels/FunctionAux.C)

Nodal AuxKernel objects abuse the notion of quadrature points, the _qp member variable is set to zero, but still must be used to access coupled variable values and material properties. This is done to allow the syntax to be consistent regardless of the AuxKernel flavor: nodal or elemental.

Mortar Nodal Auxiliary Kernel Objects

In order to compute properties in the mortar sense, it is necessary to loop over the mortar segment mesh to spatially integrate variables. MortarNodalAuxKernels offer this functionality where these "weighted" variables, which intervene in the computation of contact constraints and their residuals, can be coupled to generate the desired output value. Therefore, if postprocessing of mortar quantities is required, nodal mortar auxiliary kernels can be employed. Objects inheriting from MortarNodalAuxKernel allow for said operations on the mortar lower-dimensional domains featuring similar functionality to other nodal auxiliary kernels, including the possibility of computing quantities in an incremental manner.

Execute Flags

AuxKernel objects inherit from the SetupInterface (execute_on) so they include the "execute_on" variable. By default this parameter is set to EXEC_LINEAR and EXEC_TIMESTEP_END. The EXEC_LINEAR flag is set because it is possible to couple values compute by an AuxKernel object to other objects such as Kernel or Material objects that are used in the residual calculation. In order to ensure that the values from the auxiliary variable are correct during the iterative solve they are computed for each iteration.

However, if the auxiliary variable be computed is not being coupled to objects computing the residual evaluating the AuxKernel on each linear iteration is not necessary and can slow down the execution of a simulation. In this case, the EXEC_LINEAR flag should be removed. Likely the EXEC_INITIAL flag should be added to perform the auxiliary variable calculation during the initial setup phase as well.

Populating lower-dimensional auxiliary variables

Lower-dimensional auxiliary variables may be populated using boundary restricted auxiliary kernels. The boundary restriction of the aux kernel should be coincident with (a subset of) the lower-dimensional blocks that the lower-dimensional variable lives on. Using a boundary restricted auxiliary kernel as opposed to a lower-d block-restricted auxiliary kernel allows pulling in coincident face evaluations of higher-dimensional variables and material properties as well as evaluations of coupled lower-dimensional variables.

Example A: Post processing with AuxKernel

The following example is extracted from step 4 of the Darcy Flow and Thermomechanics Tutorial. Consider Darcy's Law for flow in porous media neglecting changes in time and gravity:

(1) where is the permeability tensor, is the fluid viscosity, and is the pressure and the velocity () may be computed as:

(2)

The left-hand side of Eq. (1) would be solved with a nonlinear variable and an appropriate Kernel object. The AuxKernel system can be used computing the velocity following Eq. (2). In the tutorial the exact calculation is performed using the DarcyVelocity object, the header and source files for this object are listed below.


#pragma once

#include "AuxKernel.h"

/**
 * Auxiliary kernel responsible for computing the Darcy velocity given
 * several fluid properties and the pressure gradient.
 */
class DarcyVelocity : public VectorAuxKernel
{
public:
  static InputParameters validParams();

  DarcyVelocity(const InputParameters & parameters);

protected:
  /**
   * AuxKernels MUST override computeValue.  computeValue() is called on
   * every quadrature point.  For Nodal Auxiliary variables those quadrature
   * points coincide with the nodes.
   */
  virtual RealVectorValue computeValue() override;

  /// The gradient of a coupled variable
  const VariableGradient & _pressure_gradient;

  /// Holds the permeability and viscosity from the material system
  const ADMaterialProperty<Real> & _permeability;
  const ADMaterialProperty<Real> & _viscosity;
};
(tutorials/darcy_thermo_mech/step04_velocity_aux/include/auxkernels/DarcyVelocity.h)

#include "DarcyVelocity.h"

#include "metaphysicl/raw_type.h"

registerMooseObject("DarcyThermoMechApp", DarcyVelocity);

InputParameters
DarcyVelocity::validParams()
{
  InputParameters params = VectorAuxKernel::validParams();

  // Add a "coupling paramater" to get a variable from the input file.
  params.addRequiredCoupledVar("pressure", "The pressure field.");

  return params;
}

DarcyVelocity::DarcyVelocity(const InputParameters & parameters)
  : VectorAuxKernel(parameters),

    // Get the gradient of the variable
    _pressure_gradient(coupledGradient("pressure")),

    // Set reference to the permeability MaterialProperty.
    // Only AuxKernels operating on Elemental Auxiliary Variables can do this
    _permeability(getADMaterialProperty<Real>("permeability")),

    // Set reference to the viscosity MaterialProperty.
    // Only AuxKernels operating on Elemental Auxiliary Variables can do this
    _viscosity(getADMaterialProperty<Real>("viscosity"))
{
}

RealVectorValue
DarcyVelocity::computeValue()
{
  // Access the gradient of the pressure at this quadrature point, then pull out the "component" of
  // it requested (x, y or z). Note, that getting a particular component of a gradient is done using
  // the parenthesis operator.
  return -MetaPhysicL::raw_value(_permeability[_qp] / _viscosity[_qp]) * _pressure_gradient[_qp];
}
(tutorials/darcy_thermo_mech/step04_velocity_aux/src/auxkernels/DarcyVelocity.C)

Example B: Decoupling Equations

Auxiliary variables may be used interchangeably with nonlinear variables with respect to coupling allowing complicated systems of equations to be decoupled for solving individually. This is very useful for testing and validation.

Consider the heat equation with an advective term that is coupled to the pressure computed in Eq. (1) as in step 6 of the Darcy Flow and Thermomechanics Tutorial:

(3) where is temperature, is the heat capacity, is the thermal conductivity, and is the porosity. The advective term () is computed in a kernel object ((tutorials/darcy_thermo_mech/step06_coupled_darcy_heat_conduction/src/kernels/DarcyAdvection.C)) and requires the pressure variable be provided as a variable:

  params.addRequiredCoupledVar("pressure", "The variable representing the pressure.");
(tutorials/darcy_thermo_mech/step06_coupled_darcy_heat_conduction/src/kernels/DarcyAdvection.C)

For testing purposes is it not desirable to include the solve for the pressure variable when examining the correctness of the heat equation solve, so an auxiliary variable that is assigned an arbitrary function of space and time is used instead. The following input file snippet demonstrates the decoupling of the pressure variable by computing it using an AuxVariable the FunctionAux object.

[AuxVariables]
  [pressure]
  []
[]

[AuxKernels]
  [pressure]
    type = FunctionAux
    variable = pressure
    function = '4000 - 3000 * x - 3000 * t*x*x*y'
    execute_on = timestep_end
  []
[]
(tutorials/darcy_thermo_mech/step06_coupled_darcy_heat_conduction/problems/step6c_decoupled.i)

Available Objects

  • Moose App
  • ADDivergenceAuxComputes the divergence of a vector of functors.
  • ADFunctorElementalAuxEvaluates a functor (variable, function or functor material property) on the current element, quadrature point, or node.
  • ADFunctorElementalGradientAuxEvaluates the gradient of a functor (variable, function or functor material property) on the current element or quadrature point.
  • ADFunctorVectorElementalAuxEvaluates a vector functor (material property usually) on the current element.For finite volume, this evaluates the vector functor at the centroid.
  • ADMaterialRankTwoTensorAuxAccess a component of a RankTwoTensor for automatic material property output
  • ADMaterialRateRealAuxOutputs element material properties rate of change
  • ADMaterialRealAuxOutputs element volume-averaged material properties
  • ADMaterialRealVectorValueAuxCapture a component of a vector material property in an auxiliary variable.
  • ADMaterialStdVectorAuxExtracts a component of a material type std::vector<Real> to an aux variable. If the std::vector is not of sufficient size then zero is returned
  • ADProjectedStatefulMaterialRankFourTensorAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ADProjectedStatefulMaterialRankTwoTensorAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ADProjectedStatefulMaterialRealAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ADProjectedStatefulMaterialRealVectorValueAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ADVectorMaterialRealVectorValueAuxConverts a vector-quantity material property into a vector auxiliary variable
  • AdvectiveFluxAuxCompute components of flux vector for advection problems .
  • ArrayParsedAuxSets field array variable values to the evaluation of a parsed expression.
  • ArrayVarReductionAuxTakes an array variable and performs a reduction operation on it (max, min, sum, average) and stores as a standard variable.
  • ArrayVariableComponentCopy a component of an array variable.
  • BuildArrayVariableAuxCombines multiple standard variables into an array variable.
  • ConstantAuxCreates a constant field in the domain.
  • ContainsPointAuxComputes a binary field where the field is 1 in the elements that contain the point and 0 everywhere else
  • CopyValueAuxReturns the specified variable as an auxiliary variable with a simple copy of the variable values.
  • DebugResidualAuxPopulate an auxiliary variable with the residual contribution of a variable.
  • DiffusionFluxAuxCompute components of flux vector for diffusion problems .
  • DivergenceAuxComputes the divergence of a vector of functors.
  • ElemExtraIDAuxPuts element extra IDs into an aux variable.
  • ElementH1ErrorFunctionAuxComputes the H1 or W^{1,p} error between an exact function and a coupled variable.
  • ElementIntegerAuxCreates a field showing the element integer.
  • ElementL2ErrorFunctionAuxA class for computing the element-wise L^2 (Euclidean) error between a function and a coupled variable.
  • ElementLengthAuxCompute the element size using Elem::hmin() or Elem::hmax() from libMesh.
  • ElementLpNormAuxCompute an elemental field variable (single value per element) equal to the Lp-norm of a coupled Variable.
  • ElementQualityAuxGenerates a field containing the quality metric for each element. Useful for visualizing mesh quality.
  • ElementUOAuxAux Kernel to display generic spatial (elemental) information from a UserObject that satisfies the underlying ElementUOProvider interface.
  • ExtraElementIDAuxPuts element extra IDs into an aux variable.
  • ForcingFunctionAuxAuxiliary Kernel that adds a forcing function to the value of an AuxVariable from the previous time step.
  • FunctionArrayAuxAuxiliary Kernel that creates and updates an array field variable by sampling functions through space and time.
  • FunctionAuxAuxiliary Kernel that creates and updates a field variable by sampling a function through space and time.
  • FunctorAuxEvaluates a functor (variable, function or functor material property) on the current element, quadrature point, or node.
  • FunctorElementalAuxEvaluates a functor (variable, function or functor material property) on the current element, quadrature point, or node.
  • FunctorElementalGradientAuxEvaluates the gradient of a functor (variable, function or functor material property) on the current element or quadrature point.
  • FunctorVectorElementalAuxEvaluates a vector functor (material property usually) on the current element.For finite volume, this evaluates the vector functor at the centroid.
  • GapValueAuxReturn the nearest value of a variable on a boundary from across a gap.
  • GhostingAuxColors the elements ghosted to the chosen PID.
  • HardwareIDAuxCreates a field showing the assignment of partitions to physical nodes in the cluster.
  • InterfaceValueUserObjectAuxGet stored value from the specified InterfaceQpUserObjectBase.
  • MaterialRankFourTensorAuxAccess a component of a RankFourTensor for automatic material property output
  • MaterialRankTwoTensorAuxAccess a component of a RankTwoTensor for automatic material property output
  • MaterialRateRealAuxOutputs element material properties rate of change
  • MaterialRealAuxOutputs element volume-averaged material properties
  • MaterialRealDenseMatrixAuxPopulate an auxiliary variable with an entry from a dense matrix material property.
  • MaterialRealTensorValueAuxObject for extracting a component of a rank two tensor material property to populate an auxiliary variable.
  • MaterialRealVectorValueAuxCapture a component of a vector material property in an auxiliary variable.
  • MaterialStdVectorAuxExtracts a component of a material type std::vector<Real> to an aux variable. If the std::vector is not of sufficient size then zero is returned
  • MaterialStdVectorRealGradientAuxExtracts a component of a material's std::vector<RealGradient> to an aux variable. If the std::vector is not of sufficient size then zero is returned
  • MeshDivisionAuxReturns the value of the mesh division index for each element / node
  • NearestNodeDistanceAuxStores the distance between a block and boundary or between two boundaries.
  • NearestNodeValueAuxRetrieves a field value from the closest node on the paired boundary and stores it on this boundary or block.
  • NodalPatchRecoveryAuxThis Auxkernel solves a least squares problem at each node to fit a value from quantities defined on quadrature points.
  • NormalizationAuxNormalizes a variable based on a Postprocessor value.
  • ParsedAuxSets a field variable value to the evaluation of a parsed expression.
  • ParsedVectorAuxSets a field vector variable value to the evaluation of a parsed expression.
  • PenetrationAuxAuxiliary Kernel for computing several geometry related quantities between two contacting bodies.
  • ProcessorIDAuxCreates a field showing the processors and partitioning.
  • ProjectedMaterialPropertyNodalPatchRecoveryAuxThis Auxkernel solves a least squares problem at each node to fit a value from quantities defined on quadrature points.
  • ProjectedStatefulMaterialRankFourTensorAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ProjectedStatefulMaterialRankTwoTensorAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ProjectedStatefulMaterialRealAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ProjectedStatefulMaterialRealVectorValueAuxPicks a component of an indexable material property to get projected on an elemental Auxvariable. For use by ProjectedStatefulMaterialStorageAction.
  • ProjectionAuxReturns the specified variable as an auxiliary variable with a projection of the source variable. If they are the same type, this amounts to a simple copy.
  • QuotientAuxDivides two coupled variables.
  • SecondTimeDerivativeAuxReturns the second order time derivative of the specified variable as an auxiliary variable.
  • SelfAuxReturns the specified variable as an auxiliary variable with a projection of the source variable. If they are the same type, this amounts to a simple copy.
  • SolutionAuxCreates fields by using information from a SolutionUserObject.
  • SpatialUserObjectAuxPopulates an auxiliary variable with a spatial value returned from a UserObject spatialValue method.
  • TagMatrixAuxCouple the diagonal of a tag matrix, and return its nodal value
  • TagVectorArrayVariableAuxCouple a tagged vector, and return its evaluations at degree of freedom indices corresponding to the coupled array variable.
  • TagVectorArrayVariableValueAuxCouple a tagged vector, and return its array value.
  • TagVectorAuxCouple a tag vector, and return its nodal value
  • TimeDerivativeAuxReturns the time derivative of the specified variable/functor as an auxiliary variable.
  • VariableGradientComponentCreates a field consisting of one component of the gradient of a coupled variable.
  • VariableTimeIntegrationAuxIntegrates a field variable in time.
  • VectorFunctionAuxAuxiliary Kernel that creates and updates a vector field variable by sampling a Function object, via the vectorValue method, through space and time.
  • VectorMagnitudeAuxCreates a field representing the magnitude of three coupled variables using an Euclidean norm.
  • VectorMaterialRealVectorValueAuxConverts a vector-quantity material property into a vector auxiliary variable
  • VectorPostprocessorVisualizationAuxRead values from a VectorPostprocessor that is producing vectors that are 'number of processors' * in length. Puts the value for each processor into an elemental auxiliary field.
  • VectorVariableComponentAuxCreates a field consisting of one component of a coupled vector variable.
  • VectorVariableMagnitudeAuxCreates a field consisting of the magnitude of a coupled vector variable.
  • VolumeAuxAuxiliary Kernel that samples volumes.
  • WeightedGapAuxReturns the specified variable as an auxiliary variable with the same value.
  • Heat Transfer App
  • JouleHeatingHeatGeneratedAuxCompute heat generated from Joule heating .
  • Peridynamics App
  • BoundaryOffsetPDClass for output offset of PD boundary nodes compared to initial FE mesh
  • NodalRankTwoPDClass for computing and outputing components and scalar quantities of nodal rank two strain and stress tensors for bond-based and ordinary state-based peridynamic models
  • NodalVolumePDClass for output nodal area(2D) or nodal volume(3D)
  • RankTwoBasedFailureCriteriaNOSPDClass for rank two tensor based failure criteria in non-ordinary state-based model
  • StretchBasedFailureCriterionPDClass for bond stretch failure criterion in bond-based model and ordinary state-based model
  • Fsi App
  • WaveHeightAuxKernelCalculates the wave heights given pressures.
  • Electromagnetics App
  • ADCurrentDensityCalculates the current density vector field (in A/m^2) when given electrostatic potential (electrostatic = true, default) or electric field.
  • CurrentDensityCalculates the current density vector field (in A/m^2) when given electrostatic potential (electrostatic = true, default) or electric field.
  • PotentialToFieldAuxAn AuxKernel that calculates the electrostatic electric field given the electrostatic potential.
  • Thermal Hydraulics App
  • ADConvectiveHeatFlux1PhaseAuxComputes convective heat flux for 1-phase flow.
  • ADVectorVelocityComponentAuxComputes the velocity from the 1D phase-fraction and area weighted momentum and density variables.
  • ConvectiveHeatFlux1PhaseAuxComputes convective heat flux for 1-phase flow.
  • MachNumberAuxComputes Mach number.
  • PrandtlNumberAuxComputes the Prandtl number for the fluid in the simulation domain
  • ReynoldsNumberAuxComputes the Reynolds number.
  • SoundSpeedAuxComputes the speed of sound.
  • SpecificTotalEnthalpyAuxCompute the specific total enthalpy
  • SumAuxSum of nonlinear or auxiliary variables
  • THMSpecificInternalEnergyAuxComputed the specific internal energy.
  • THMSpecificVolumeAuxComputes the specific volume for the phase.
  • VariableValueTransferAuxRetrieves a field value from the closest node on the paired boundary and stores it on this boundary or block.
  • VectorVelocityComponentAuxComputes the component of a vector-valued velocity field given by its magnitude and direction.
  • WeightedAverageAuxWeighted average of variables using other variables as weights
  • Contact App
  • CohesiveZoneMortarUserObjectAuxPopulates an auxiliary variable with mortar cohesive zone model quantities.
  • ContactPressureAuxComputes the contact pressure from the contact force and nodal area
  • MortarArchardsLawAuxReturns the weighted gap velocity at a node. This quantity is useful for mortar contact, particularly when dual basis functions are used in contact mechanics
  • MortarFrictionalPressureVectorAuxThis class creates an auxiliary vector for outputting the mortar frictional pressure vector.
  • MortarFrictionalStateAuxThis class creates discrete states for nodes into frictional contact, including contact/no-contact and stick/slip.
  • MortarPressureComponentAuxThis class transforms the Cartesian Lagrange multiplier vector to local coordinates and outputs each individual component along the normal or tangential direction.
  • PenaltyMortarUserObjectAuxPopulates an auxiliary variable with a contact quantities from penalty mortar contact.
  • WeightedGapVelAuxReturns the weighted gap velocity at a node. This quantity is useful for mortar contact, particularly when dual basis functions are used in contact mechanics
  • XFEMApp
  • CutSubdomainIDAuxFill the elemental variable with CutSubdomainID
  • MeshCutLevelSetAuxCalculates signed distance from interface defined by InterfaceMeshCutUserObject.
  • XFEMCutPlaneAuxComputes the normal and origin of a cutting plane for each partial element.
  • XFEMMarkerAuxIdentify the crack tip elements.
  • XFEMVolFracAuxComputes the volume fraction of the physical material in each partial element.
  • Solid Mechanics App
  • ADKineticEnergyAuxCompute the kinetic energy of continuum-based finite elements
  • ADRankFourAuxAccess a component of a RankFourTensor
  • ADRankTwoAuxAccess a component of a RankTwoTensor
  • ADRankTwoScalarAuxCompute a scalar property of a RankTwoTensor
  • AccumulateAux
  • CylindricalRankTwoAuxTakes RankTwoTensor material and outputs component in cylindrical coordinates
  • DomainIntegralQFunctionComputes the q-function for a segment along the crack front, used in the calculation of the J-integral
  • DomainIntegralTopologicalQFunctionDetermines if a node is within the ring of the crack front defintion; this object is normally created by the DomainIntegralAction.
  • ElasticEnergyAuxCompute the local elastic energy
  • GlobalDisplacementAuxAuxKernel to visualize the displacements generated by the global strain tensor
  • KineticEnergyAuxCompute the kinetic energy of continuum-based finite elements
  • NewmarkAccelAuxComputes the current acceleration using the Newmark method.
  • NewmarkVelAuxCalculates the current velocity using Newmark method.
  • RadialDisplacementCylinderAuxCompute the radial component of the displacement vector for cylindrical models.
  • RadialDisplacementSphereAuxCompute the radial component of the displacement vector for spherical models.
  • RankFourAuxAccess a component of a RankFourTensor
  • RankTwoAuxAccess a component of a RankTwoTensor
  • RankTwoScalarAuxCompute a scalar property of a RankTwoTensor
  • RotationAngleCompute the field of angular rotations of points around an axis defined by an origin point and a direction vector
  • TestNewmarkTIAssigns the velocity/acceleration calculated by time integrator to the velocity/acceleration auxvariable.
  • Phase Field App
  • BndsCalcAuxCalculate location of grain boundaries in a polycrystalline sample
  • CrossTermGradientFreeEnergyFree energy contribution from the cross terms in ACMultiInterface
  • DiscreteNucleationAuxProject the DiscreteNucleationMap state onto an AuxVariable
  • EBSDReaderAvgDataAux
  • EBSDReaderPointDataAux
  • EulerAngleProvider2RGBAuxOutput RGB representation of crystal orientation from user object to an AuxVariable. The entire domain must have the same crystal structure.
  • EulerAngleVariables2RGBAux
  • FeatureFloodCountAuxFeature detection by connectivity analysis
  • GrainAdvectionAuxCalculates the advection velocity of grain due to rigid body translation and rotation
  • GrainBoundaryVelocityCompute the velocity of grain boundaries.
  • KKSGlobalFreeEnergyTotal free energy in KKS system, including chemical, barrier and gradient terms
  • KKSMultiFreeEnergyTotal free energy in multi-phase KKS system, including chemical, barrier and gradient terms
  • LinearizedInterfaceAuxCalculates the order parameter from the linearized interface function
  • OutputEulerAnglesOutput Euler angles from user object to an AuxVariable.
  • PFCEnergyDensity
  • PFCRFFEnergyDensity
  • SolutionAuxMisorientationBoundaryCalculate location of grain boundaries by using information from a SolutionUserObject.
  • TotalFreeEnergyTotal free energy (both the bulk and gradient parts), where the bulk free energy has been defined in a material
  • Stochastic Tools App
  • SurrogateModelArrayAuxKernelSets a value of a variable based on a surrogate model.
  • SurrogateModelAuxKernelSets a value of a variable based on a surrogate model.
  • Navier Stokes App
  • CourantComputes |u| dt / h_min.
  • EnthalpyAuxThis AuxKernel computes the specific enthalpy of the fluidfrom the total energy and the pressure.
  • HasPorosityJumpFaceShows whether an element has any attached porosity jump faces
  • INSCourantComputes h_min / |u|.
  • INSFVMixingLengthTurbulentViscosityAuxComputes the turbulent viscosity for the mixing length model.
  • INSQCriterionAuxThis class computes the Q criterion, a scalar whichaids in vortex identification in turbulent flows
  • INSStressComponentAuxThis class computes the stress component based on pressure and velocity for incompressible Navier-Stokes
  • InternalEnergyAuxThis AuxKernel computes the internal energy based on the equation of state / fluid properties and the local pressure and density.
  • NSInternalEnergyAuxAuxiliary kernel for computing the internal energy of the fluid.
  • NSLiquidFractionAuxComputes liquid fraction given the temperature.
  • NSMachAuxAuxiliary kernel for computing the Mach number assuming an ideal gas.
  • NSPressureAuxNodal auxiliary variable, for computing pressure at the nodes.
  • NSSpecificTotalEnthalpyAuxNodal auxiliary variable, for computing enthalpy at the nodes.
  • NSTemperatureAuxTemperature is an auxiliary value computed from the total energy based on the FluidProperties.
  • NSVelocityAuxVelocity auxiliary value.
  • PecletNumberFunctorAuxComputes the Peclet number: u*L/alpha.
  • ReynoldsNumberFunctorAuxComputes rho*u*L/mu.
  • SpecificInternalEnergyAuxThis AuxKernel computes the specific internal energy based from the total and the kinetic energy.
  • SpecificVolumeAuxThis auxkernel computes the specific volume of the fluid.
  • TurbulentConductivityAuxCalculates the turbulent heat conductivity.
  • WallDistanceMixingLengthAuxComputes the turbulent mixing length by assuming that it is proportional to the distance from the nearest wall. The mixinglength is capped at a distance proportional to inputted parameter delta.
  • WallFunctionWallShearStressAuxCalculates the wall shear stress based on algebraic standard velocity wall functions.
  • WallFunctionYPlusAuxCalculates y+ value according to the algebraic velocity standard wall function.
  • kEpsilonViscosityAuxCalculates the turbulent viscosity according to the k-epsilon model.
  • Misc App
  • CoupledDirectionalMeshHeightInterpolationScales a variable based on position relative to the model bounds in a specified direction
  • Geochemistry App
  • GeochemistryQuantityAuxExtracts information from the Reactor and records it in the AuxVariable
  • NodalVoidVolumeAuxExtracts information from the NodalVoidVolume UserObject and records it in the AuxVariable
  • Porous Flow App
  • ADPorousFlowDarcyVelocityComponentDarcy velocity (in m3.s-1.m-2, or m.s-1) -(k_ij * krel /mu (nabla_j P - w_j)), where k_ij is the permeability tensor, krel is the relative permeability, mu is the fluid viscosity, P is the fluid pressure, and w_j is the fluid weight.
  • ADPorousFlowPropertyAuxAuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.
  • PorousFlowDarcyVelocityComponentDarcy velocity (in m3.s-1.m-2, or m.s-1) -(k_ij * krel /mu (nabla_j P - w_j)), where k_ij is the permeability tensor, krel is the relative permeability, mu is the fluid viscosity, P is the fluid pressure, and w_j is the fluid weight.
  • PorousFlowDarcyVelocityComponentLowerDimensionalDarcy velocity on a lower-dimensional element embedded in a higher-dimensional mesh. Units m3.s-1.m-2, or m.s-1. Darcy velocity = -(k_ij * krel /(mu * a) (nabla_j P - w_j)), where k_ij is the permeability tensor, krel is the relative permeability, mu is the fluid viscosity, P is the fluid pressure, a is the fracture aperture and w_j is the fluid weight. The difference between this AuxKernel and PorousFlowDarcyVelocity is that this one projects gravity along the element's tangent direction. NOTE! For a meaningful answer, your permeability tensor must NOT contain terms that rotate tangential vectors to non-tangential vectors.
  • PorousFlowElementLengthAuxKernel to compute the 'length' of elements along a given direction. A plane is constructed through the element's centroid, with normal equal to the direction given. The average of the distance of the nodal positions to this plane is the 'length' returned. The Variable for this AuxKernel must be an elemental Variable
  • PorousFlowElementNormalAuxKernel to compute components of the element normal. This is mostly designed for 2D elements living in 3D space, however, the 1D-element and 3D-element cases are handled as special cases. The Variable for this AuxKernel must be an elemental Variable
  • PorousFlowPropertyAuxAuxKernel to provide access to properties evaluated at quadpoints. Note that elemental AuxVariables must be used, so that these properties are integrated over each element.
  • Functional Expansion Tools App
  • FunctionSeriesToAuxAuxKernel to convert a functional expansion (Functions object, type = FunctionSeries) to an AuxVariable
  • Fluid Properties App
  • FluidDensityAuxComputes density from pressure and temperature
  • PressureAuxComputes pressure given specific volume and specific internal energy
  • SaturationTemperatureAuxComputes saturation temperature from pressure and 2-phase fluid properties object
  • SpecificEnthalpyAuxComputes specific enthalpy from pressure and temperature
  • StagnationPressureAuxComputes stagnation pressure from specific volume, specific internal energy, and velocity
  • StagnationTemperatureAuxComputes stagnation temperature from specific volume, specific internal energy, and velocity
  • TemperatureAuxComputes temperature given specific volume and specific internal energy
  • Chemical Reactions App
  • AqueousEquilibriumRxnAuxConcentration of secondary equilibrium species
  • EquilibriumConstantAuxEquilibrium constant for a given equilibrium species (in form log10(Keq))
  • KineticDisPreConcAuxConcentration of secondary kinetic species
  • KineticDisPreRateAuxKinetic rate of secondary kinetic species
  • PHAuxpH of solution
  • TotalConcentrationAuxTotal concentration of primary species (including stoichiometric contribution to secondary equilibrium species)

Available Subsystems

Available Actions