Kokkos BCs System
Before reading this documentation, consider reading the following materials first for a better understanding of this documentation:
BCs System to understand the MOOSE boundary condition system,
Getting Started with Kokkos-MOOSE to understand the programming practices for Kokkos-MOOSE,
Kokkos Kernels System to understand the common design pattern of objects in Kokkos-MOOSE.
Kokkos-MOOSE boundary conditions do not support automatic differention yet.
The basic design pattern of Kokkos-MOOSE kernels described in Kokkos Kernels System applies to the boundary conditions as well. You can create your own integrated and nodal boundary conditions by subclassing Moose::Kokkos::IntegratedBC and Moose::Kokkos::NodalBC, respectively, and following the same pattern with kernels including registering your boundary conditions with registerKokkosResidualObject(). Especially, integrated boundary conditions have identical interfaces with kernels, so they will not be explained here in detail. See the following source codes of KokkosCoupledVarNeumannBC for an example of an integrated boundary condition:
Listing 1: The KokkosCoupledVarNeumannBC header file.
#pragma once
#include "KokkosIntegratedBC.h"
/**
* Implements a Neumann BC where grad(u)=_coupled_var on the boundary.
* Uses the term produced from integrating the diffusion operator by parts.
*/
class KokkosCoupledVarNeumannBC : public Moose::Kokkos::IntegratedBC
{
public:
static InputParameters validParams();
KokkosCoupledVarNeumannBC(const InputParameters & parameters);
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
const unsigned int qp,
AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int i,
const unsigned int j,
const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const;
protected:
/// Variable providing the value of grad(u) on the boundary.
const Moose::Kokkos::VariableValue _coupled_var;
/// The identifying number of the coupled variable
const unsigned int _coupled_num;
/// A coefficient that is multiplied with the residual contribution
const Real _coef;
/// Scale factor
const Moose::Kokkos::VariableValue _scale_factor;
};
KOKKOS_FUNCTION inline Real
KokkosCoupledVarNeumannBC::computeQpResidual(const unsigned int i,
const unsigned int qp,
AssemblyDatum & datum) const
{
return -_scale_factor(datum, qp) * _coef * _test(datum, i, qp) * _coupled_var(datum, qp);
}
KOKKOS_FUNCTION inline Real
KokkosCoupledVarNeumannBC::computeQpOffDiagJacobian(const unsigned int i,
const unsigned int j,
const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const
{
if (jvar == _coupled_num)
return -_scale_factor(datum, qp) * _coef * _test(datum, i, qp) * _phi(datum, j, qp);
else
return 0;
}
(framework/include/kokkos/bcs/KokkosCoupledVarNeumannBC.h)Listing 2: The KokkosCoupledVarNeumannBC source file.
#include "KokkosCoupledVarNeumannBC.h"
registerKokkosResidualObject("MooseApp", KokkosCoupledVarNeumannBC);
InputParameters
KokkosCoupledVarNeumannBC::validParams()
{
InputParameters params = IntegratedBC::validParams();
params.addRequiredCoupledVar("v", "Coupled variable setting the gradient on the boundary.");
params.addCoupledVar("scale_factor", 1., "Scale factor to multiply the heat flux with");
params.addParam<Real>(
"coef", 1.0, "Coefficent ($\\sigma$) multiplier for the coupled force term.");
params.addClassDescription("Imposes the integrated boundary condition "
"$\\frac{\\partial u}{\\partial n}=v$, "
"where $v$ is a variable.");
return params;
}
KokkosCoupledVarNeumannBC::KokkosCoupledVarNeumannBC(const InputParameters & parameters)
: IntegratedBC(parameters),
_coupled_var(kokkosCoupledValue("v")),
_coupled_num(coupled("v")),
_coef(getParam<Real>("coef")),
_scale_factor(kokkosCoupledValue("scale_factor"))
{
}
(framework/src/kokkos/bcs/KokkosCoupledVarNeumannBC.K)On the other hand, nodal boundary conditions have slightly different interfaces. The hook methods for a nodal boundary condition have the following signatures:
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int qp, AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const;
The test and trial function indices, i and j, are no longer passed as arguments. qp corresponds to the dummy _qp index in the original MOOSE and is always zero. To keep the consistency between interfaces, however, it is still passed as an argument and used for getting solution values. _current_node, which is a pointer to the current libMesh node object, does not have a direct replacement. Instead, the node index can be queried by datum.node() and used to retrieve mesh data from the Kokkos mesh object. The node coordinate can also be obtained by datum.q_point(qp). As a result, the following residual function in DirichletBCBase:
Real
DirichletBCBase::computeQpResidual()
{
return _u[_qp] - computeQpValue();
}
becomes the following in Moose::Kokkos::DirichletBCBase:
template <typename Derived>
KOKKOS_FUNCTION Real
DirichletBCBase<Derived>::computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const
{
auto bc = static_cast<const Derived *>(this);
return _u(datum, qp) - bc->computeValue(qp, datum);
}
Also note here the static implementation of computeValue using the Curiously Recurring Template Pattern (CRTP) which is originally a virtual function.
See the following source codes of KokkosMatchedValueBC for another example of a nodal boundary condition:
Listing 3: The KokkosMatchedValueBC header file.
#pragma once
#include "KokkosNodalBC.h"
/**
* Implements a simple coupled boundary condition where u=v on the boundary.
*/
class KokkosMatchedValueBC : public Moose::Kokkos::NodalBC
{
public:
static InputParameters validParams();
KokkosMatchedValueBC(const InputParameters & parameters);
KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const;
KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int jvar,
const unsigned int qp,
AssemblyDatum & datum) const;
protected:
const Moose::Kokkos::VariableValue _v;
/// The id of the coupled variable
const unsigned int _v_num;
/// Coefficient for primary variable
const Real _u_coeff;
/// Coefficient for coupled variable
const Real _v_coeff;
};
KOKKOS_FUNCTION inline Real
KokkosMatchedValueBC::computeQpResidual(const unsigned int qp, AssemblyDatum & datum) const
{
return _u_coeff * _u(datum, qp) - _v_coeff * _v(datum, qp);
}
KOKKOS_FUNCTION inline Real
KokkosMatchedValueBC::computeQpOffDiagJacobian(const unsigned int jvar,
const unsigned int /* qp */,
AssemblyDatum & /* datum */) const
{
if (jvar == _v_num)
return -_v_coeff;
else
return 0;
}
(framework/include/kokkos/bcs/KokkosMatchedValueBC.h)Listing 4: The KokkosMatchedValueBC source file.
#include "KokkosMatchedValueBC.h"
registerKokkosResidualObject("MooseApp", KokkosMatchedValueBC);
InputParameters
KokkosMatchedValueBC::validParams()
{
InputParameters params = NodalBC::validParams();
params.addRequiredCoupledVar("v", "The variable whose value we are to match.");
params.addParam<Real>("u_coeff", 1.0, " A coefficient for primary variable u");
params.addParam<Real>("v_coeff", 1.0, " A coefficient for coupled variable v");
params.addClassDescription("Implements a NodalBC which equates two different Variables' values "
"on a specified boundary.");
return params;
}
KokkosMatchedValueBC::KokkosMatchedValueBC(const InputParameters & parameters)
: NodalBC(parameters),
_v(kokkosCoupledNodalValue("v")),
_v_num(coupled("v")),
_u_coeff(getParam<Real>("u_coeff")),
_v_coeff(getParam<Real>("v_coeff"))
{
}
(framework/src/kokkos/bcs/KokkosMatchedValueBC.K)Available Objects
- Moose App
- ADConservativeAdvectionBCBoundary condition for advection when it is integrated by parts. Supports Dirichlet (inlet-like) and implicit (outlet-like) conditions.
- ADCoupledVarNeumannBCImposes the integrated boundary condition , where is a variable.
- ADDirichletBCImposes the essential boundary condition , where is a constant, controllable value.
- ADFunctionDirichletBCImposes the essential boundary condition , where is calculated by a function.
- ADFunctionNeumannBCImposes the integrated boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- ADFunctionPenaltyDirichletBCEnforces a (possibly) time and space-dependent MOOSE Function Dirichlet boundary condition in a weak sense by penalizing differences between the current solution and the Dirichlet data.
- ADMatNeumannBCImposes the integrated boundary condition , where is a constant, is a material property, and is a coefficient defined by the kernel for .
- ADMatchedValueBCImplements a NodalBC which equates two different Variables' values on a specified boundary.
- ADNeumannBCImposes the integrated boundary condition , where is a constant, controllable value.
- ADPenaltyDirichletBCEnforces a Dirichlet boundary condition in a weak sense by penalizing differences between the current solution and the Dirichlet data.
- ADRobinBCImposes the Robin integrated boundary condition .
- ADVectorFunctionDirichletBCImposes the essential boundary condition , where components are calculated with functions.
- ADVectorFunctionNeumannBCImposes the integrated boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- ADVectorMatchedValueBCImplements a ADVectorNodalBC which equates two different Variables' values on a specified boundary.
- ADVectorRobinBCImposes the Robin integrated boundary condition .
- AdvectionIPHDGDirichletBCWeakly imposes Dirichlet boundary conditions for a hybridized discretization of an advection equation
- AdvectionIPHDGOutflowBCImplements an outflow boundary condition for use with a hybridized discretization of the advection equation
- AdvectionIPHDGPrescribedFluxBCImplements a prescribed flux condition for use with a hybridized discretization of the advection equation
- ArrayDirichletBCImposes the essential boundary condition , where are constant, controllable values.
- ArrayHFEMDirichletBCImposes the Dirichlet BC with HFEM.
- ArrayNeumannBCImposes the integrated boundary condition , where is a constant, controllable value.
- ArrayPenaltyDirichletBCEnforces a Dirichlet boundary condition in a weak sense with , where is the constant scalar penalty; is the test functions and is the differences between the current solution and the Dirichlet data.
- ArrayVacuumBCImposes the Robin boundary condition .
- ConservativeAdvectionBCBoundary condition for advection when it is integrated by parts. Supports Dirichlet (inlet-like) and implicit (outlet-like) conditions.
- ConvectiveFluxBCDetermines boundary values via the initial and final values, flux, and exposure duration
- CoupledVarNeumannBCImposes the integrated boundary condition , where is a variable.
- DGFunctionDiffusionDirichletBCDiffusion Dirichlet boundary condition for discontinuous Galerkin method.
- DiffusionFluxBCComputes a boundary residual contribution consistent with the Diffusion Kernel. Does not impose a boundary condition; instead computes the boundary contribution corresponding to the current value of grad(u) and accumulates it in the residual vector.
- DiffusionIPHDGDirichletBCWeakly imposes Dirichlet boundary conditions for a hybridized discretization of a diffusion equation
- DiffusionIPHDGPrescribedFluxBCImplements a flux boundary condition for use with a hybridized discretization of the diffusion equation
- DiffusionLHDGDirichletBCWeakly imposes Dirichlet boundary conditions for a hybridized discretization of a diffusion equation
- DiffusionLHDGPrescribedGradientBCImplements a flux boundary condition for use with a hybridized discretization of the diffusion equation
- DirectionalNeumannBCImposes the integrated boundary condition , where is a user-defined, constant vector.
- DirichletBCImposes the essential boundary condition , where is a constant, controllable value.
- EigenArrayDirichletBCArray Dirichlet BC for eigenvalue solvers
- EigenDirichletBCDirichlet BC for eigenvalue solvers
- FunctionDirichletBCImposes the essential boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- FunctionGradientNeumannBCImposes the integrated boundary condition arising from integration by parts of a diffusion/heat conduction operator, and where the exact solution can be specified.
- FunctionNeumannBCImposes the integrated boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- FunctionPenaltyDirichletBCEnforces a (possibly) time and space-dependent MOOSE Function Dirichlet boundary condition in a weak sense by penalizing differences between the current solution and the Dirichlet data.
- FunctorDirichletBCImposes the Dirichlet boundary condition , where is a functor and can have complex dependencies.
- FunctorNeumannBCImposes the integrated boundary condition , where is a functor.
- HFEMDirichletBCImposes the Dirichlet BC with HFEM.
- KokkosCoupledVarNeumannBCImposes the integrated boundary condition , where is a variable.
- KokkosDirichletBCImposes the essential boundary condition , where is a constant, controllable value.
- KokkosMatchedValueBCImplements a NodalBC which equates two different Variables' values on a specified boundary.
- KokkosNeumannBCImposes the integrated boundary condition , where is a constant, controllable value.
- KokkosPostprocessorNeumannBCNeumann boundary condition with value prescribed by a Postprocessor value.
- KokkosVacuumBCVacuum boundary condition for diffusion.
- LagrangeVecDirichletBCImposes the essential boundary condition , where are constant, controllable values.
- LagrangeVecFunctionDirichletBCImposes the essential boundary condition , where components are calculated with functions.
- MFEMBoundaryIntegratedBCAdds the boundary integrator to an MFEM problem for the linear form arising from the weak form of the forcing term .
- MFEMBoundaryNormalIntegratedBCAdds the boundary integrator to an MFEM problem for the linear form
- MFEMConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties to add to MFEM problems.
- MFEMScalarDirichletBCApplies a Dirichlet condition to a scalar variable.
- MFEMVectorBoundaryIntegratedBCAdds the boundary integrator to an MFEM problem for the linear form
- MFEMVectorDirichletBCApplies a Dirichlet condition to all components of a vector variable.
- MFEMVectorFEBoundaryFluxIntegratedBCAdds the boundary integrator to an MFEM problem for the linear form
- MFEMVectorNormalDirichletBCApplies a Dirichlet condition to the normal components of a vector variable.
- MFEMVectorTangentialDirichletBCApplies a Dirichlet condition to the tangential components of a vector variable.
- MatNeumannBCImposes the integrated boundary condition , where is a constant, is a material property, and is a coefficient defined by the kernel for .
- MatchedValueBCImplements a NodalBC which equates two different Variables' values on a specified boundary.
- NeumannBCImposes the integrated boundary condition , where is a constant, controllable value.
- OneDEqualValueConstraintBCComputes the integral of lambda times dg term from the mortar method (for two 1D domains only).
- PenaltyDirichletBCEnforces a Dirichlet boundary condition in a weak sense by penalizing differences between the current solution and the Dirichlet data.
- PostprocessorDirichletBCDirichlet boundary condition with value prescribed by a Postprocessor value.
- PostprocessorNeumannBCNeumann boundary condition with value prescribed by a Postprocessor value.
- SinDirichletBCImposes a time-varying essential boundary condition , where varies from an given initial value at time to a given final value over a specified duration.
- SinNeumannBCImposes a time-varying flux boundary condition , where varies from an given initial value at time to a given final value over a specified duration.
- VacuumBCVacuum boundary condition for diffusion.
- VectorCurlPenaltyDirichletBCEnforces a Dirichlet boundary condition for the curl of vector nonlinear variables in a weak sense by applying a penalty to the difference in the current solution and the Dirichlet data.
- VectorDirichletBCImposes the essential boundary condition , where are constant, controllable values.
- VectorDivPenaltyDirichletBCEnforces, in a weak sense, a Dirichlet boundary condition on the divergence of a nonlinear vector variable by applying a penalty to the difference between the current solution and the Dirichlet data.
- VectorFunctionDirichletBCImposes the essential boundary condition , where components are calculated with functions.
- VectorNeumannBCImposes the integrated boundary condition , where is a user-defined, constant vector.
- VectorPenaltyDirichletBCEnforces a Dirichlet boundary condition for vector nonlinear variables in a weak sense by applying a penalty to the difference in the current solution and the Dirichlet data.
- WeakGradientBCComputes a boundary residual contribution consistent with the Diffusion Kernel. Does not impose a boundary condition; instead computes the boundary contribution corresponding to the current value of grad(u) and accumulates it in the residual vector.
- XFEMApp
- CrackTipEnrichmentCutOffBCImposes the essential boundary condition , where is a constant, controllable value.
- Scalar Transport App
- BinaryRecombinationBCModels recombination of the variable with a coupled species at boundaries, resulting in loss
- DissociationFluxBCModels creation of the variable at boundaries due to dissociation of a coupled variable, e.g. B -> A
- Porous Flow App
- PorousFlowEnthalpySinkApplies a source equal to the product of the mass flux and the fluid enthalpy. The enthalpy is computed at temperature T_in and pressure equal to the porepressure in the porous medium, if fluid_phase is given, otherwise at the supplied porepressure. Hence this adds heat energy to the porous medium at rate corresponding to a fluid being injected at (porepressure, T_in) at rate (-flux_function).
- PorousFlowHalfCubicSinkApplies a flux sink to a boundary. The base flux defined by PorousFlowSink is multiplied by a cubic.
- PorousFlowHalfGaussianSinkApplies a flux sink to a boundary. The base flux defined by PorousFlowSink is multiplied by a Gaussian.
- PorousFlowOutflowBCApplies an 'outflow' boundary condition, which allows fluid components or heat energy to flow freely out of the boundary as if it weren't there. This is fully upwinded
- PorousFlowPiecewiseLinearSinkApplies a flux sink to a boundary. The base flux defined by PorousFlowSink is multiplied by a piecewise linear function.
- PorousFlowSinkApplies a flux sink to a boundary.
- Functional Expansion Tools App
- FXFluxBCSets a flux boundary condition, evaluated using a FunctionSeries instance. This does not fix the flux, but rather 'strongly encourages' flux agreement by penalizing the differences through contributions to the residual.
- FXValueBCImposes a fixed value boundary condition, evaluated using a FunctionSeries instance.
- FXValuePenaltyBCSets a value boundary condition, evaluated using a FunctionSeries instance. This does not fix the value, but rather 'strongly encourages' value agreement by penalizing the differences through contributions to the residual.
- Navier Stokes App
- AdvectionBCBoundary conditions for outflow/outflow of advected quantities: phi * velocity * normal, where phi is the advected quantitiy
- EnergyFreeBCImplements free advective flow boundary conditions for the energy equation.
- FluidWallMomentumBCImplicitly sets normal component of velocity to zero if the advection term of the momentum equation is integrated by parts
- INSADDisplaceBoundaryBCBoundary condition for displacing a boundary
- INSADDummyDisplaceBoundaryIntegratedBCThis object adds Jacobian entries for the boundary displacement dependence on the velocity
- INSADMomentumNoBCBCThis class implements the 'No BC' boundary condition based on the 'Laplace' form of the viscous stress tensor.
- INSADSurfaceTensionBCSurface tension stresses.
- INSADVaporRecoilPressureMomentumFluxBCVapor recoil pressure momentum flux
- INSFEFluidEnergyBCSpecifies flow of energy through a boundary
- INSFEFluidEnergyDirichletBCImposes a Dirichlet condition on temperature at inlets. Is not applied at outlets
- INSFEFluidMassBCSpecifies flow of mass through a boundary given a velocity function or postprocessor
- INSFEFluidMomentumBCSpecifies flow of momentum through a boundary
- INSFEFluidWallMomentumBCImplicitly sets normal component of velocity to zero if the advection term of the momentum equation is integrated by parts
- INSFEMomentumFreeSlipBCImplements free slip boundary conditions for the Navier Stokesmomentum equation.
- INSMomentumNoBCBCLaplaceFormThis class implements the 'No BC' boundary condition based on the 'Laplace' form of the viscous stress tensor.
- INSMomentumNoBCBCTractionFormThis class implements the 'No BC' boundary condition based on the 'traction' form of the viscous stress tensor.
- INSTemperatureNoBCBCThis class implements the 'No BC' boundary condition discussed by Griffiths, Papanastiou, and others.
- ImplicitNeumannBCThis class implements a form of the Neumann boundary condition in which the boundary term is treated 'implicitly'.
- MDFluidEnergyBCSpecifies flow of energy through a boundary
- MDFluidEnergyDirichletBCImposes a Dirichlet condition on temperature at inlets. Is not applied at outlets
- MDFluidMassBCSpecifies flow of mass through a boundary given a velocity function or postprocessor
- MDFluidMomentumBCSpecifies flow of momentum through a boundary
- MDMomentumFreeSlipBCImplements free slip boundary conditions for the Navier Stokesmomentum equation.
- MassFluxPenaltyBCAdds the exterior boundary contribution of penalized jumps in the velocity variable in one component of the momentum equations.
- MassFreeBCImplements free advective flow boundary conditions for the mass equation.
- MassMatrixIntegratedBCComputes a finite element mass matrix meant for use in preconditioning schemes which require one
- MomentumFreeBCImplements free flow boundary conditions for one of the momentum equations.
- MomentumFreeSlipBCImplements free slip boundary conditions for the Navier Stokesmomentum equation.
- NSEnergyInviscidSpecifiedBCThis class corresponds to the inviscid part of the 'natural' boundary condition for the energy equation.
- NSEnergyInviscidSpecifiedDensityAndVelocityBCThis class corresponds to the inviscid part of the 'natural' boundary condition for the energy equation.
- NSEnergyInviscidSpecifiedNormalFlowBCThis class corresponds to the inviscid part of the 'natural' boundary condition for the energy equation.
- NSEnergyInviscidSpecifiedPressureBCThis class corresponds to the inviscid part of the 'natural' boundary condition for the energy equation.
- NSEnergyInviscidUnspecifiedBCThis class corresponds to the inviscid part of the 'natural' boundary condition for the energy equation.
- NSEnergyViscousBCThis class couples together all the variables for the compressible Navier-Stokes equations to allow them to be used in derived IntegratedBC classes.
- NSEnergyWeakStagnationBCThe inviscid energy BC term with specified normal flow.
- NSImposedVelocityBCImpose Velocity BC.
- NSImposedVelocityDirectionBCThis class imposes a velocity direction component as a Dirichlet condition on the appropriate momentum equation.
- NSInflowThermalBCThis class is used on a boundary where the incoming flow values (rho, u, v, T) are all completely specified.
- NSMassSpecifiedNormalFlowBCThis class implements the mass equation boundary term with a specified value of rho*(u.n) imposed weakly.
- NSMassUnspecifiedNormalFlowBCThis class implements the mass equation boundary term with the rho*(u.n) boundary integral computed implicitly.
- NSMassWeakStagnationBCThe inviscid energy BC term with specified normal flow.
- NSMomentumConvectiveWeakStagnationBCThe convective part (sans pressure term) of the momentum equation boundary integral evaluated at specified stagnation temperature, stagnation pressure, and flow direction values.
- NSMomentumInviscidNoPressureImplicitFlowBCMomentum equation boundary condition used when pressure is not integrated by parts.
- NSMomentumInviscidSpecifiedNormalFlowBCMomentum equation boundary condition in which pressure is specified (given) and the value of the convective part is allowed to vary (is computed implicitly).
- NSMomentumInviscidSpecifiedPressureBCMomentum equation boundary condition in which pressure is specified (given) and the value of the convective part is allowed to vary (is computed implicitly).
- NSMomentumPressureWeakStagnationBCThis class implements the pressure term of the momentum equation boundary integral for use in weak stagnation boundary conditions.
- NSMomentumViscousBCThis class corresponds to the viscous part of the 'natural' boundary condition for the momentum equations.
- NSPenalizedNormalFlowBCThis class penalizes the the value of u.n on the boundary so that it matches some desired value.
- NSPressureNeumannBCThis kernel is appropriate for use with a 'zero normal flow' boundary condition in the context of the Euler equations.
- NSStagnationPressureBCThis Dirichlet condition imposes the condition p_0 = p_0_desired.
- NSStagnationTemperatureBCThis Dirichlet condition imposes the condition T_0 = T_0_desired.
- NSThermalBCNS thermal BC.
- NavierStokesLHDGOutflowBCImplements an outflow boundary condition for use with a hybridized discretization of the incompressible Navier-Stokes equations
- NavierStokesLHDGVelocityDirichletBCWeakly imposes Dirichlet boundary conditions for the velocity for a hybridized discretization of the Navier-Stokes equations
- NavierStokesStressIPHDGDirichletBCWeakly imposes Dirichlet boundary conditions for a hybridized discretization of a Navier-Stokes equation stress term
- Heat Transfer App
- ADConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficient given by material properties.
- ADFunctionRadiativeBCBoundary condition for radiative heat exchange where the emissivity function is supplied by a Function.
- ADInfiniteCylinderRadiativeBCBoundary condition for radiative heat exchange with a cylinderwhere the boundary is approximated as a cylinder as well.
- ConvectiveFluxFunctionDetermines boundary value by fluid heat transfer coefficient and far-field temperature
- ConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties.
- CoupledConvectiveFluxIntegrated boundary condition for modeling a convective heat flux.
- CoupledConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficent given by auxiliary variables.
- DirectionalFluxBCApplies a directional flux multiplied by the surface normal vector. Can utilize the self shadowing calculation from a SelfShadowSideUserObject.
- FunctionRadiativeBCBoundary condition for radiative heat exchange where the emissivity function is supplied by a Function.
- GapHeatTransferTransfers heat across a gap between two surfaces dependent on the gap geometry specified.
- GapPerfectConductanceEnforces equal temperatures across the gap.
- GaussianEnergyFluxBCDescribes an incoming heat flux beam with a Gaussian profile
- GrayLambertNeumannBCThis BC imposes a heat flux density that is computed from the GrayLambertSurfaceRadiationBase userobject.
- HeatConductionBCBoundary condition for a diffusive heat flux.
- InfiniteCylinderRadiativeBCBoundary condition for radiative heat exchange with a cylinderwhere the boundary is approximated as a cylinder as well.
- KokkosConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficent given by material properties.
- KokkosCoupledConvectiveHeatFluxBCConvective heat transfer boundary condition with temperature and heat transfer coefficent given by auxiliary variables.
- Rdg App
- AEFVBCA boundary condition kernel for the advection equation using a cell-centered finite volume method.
- Peridynamics App
- RBMPresetOldValuePDClass to apply a preset BC to nodes with rigid body motion (RBM).
- Electromagnetics App
- EMRobinBCFirst order Robin-style Absorbing/Port BC for scalar variables, assuming plane waves.
- VectorEMRobinBCFirst order Robin-style Absorbing/Port BC for vector variables.
- VectorTransientAbsorbingBCFirst order transient absorbing boundary condition for vector variables.
- Fsi App
- FluidFreeSurfaceBCApplies a mixed Dirichlet-Neumann BC on the fluid surface.
- Solid Mechanics App
- ADPenaltyInclinedNoDisplacementBCPenalty Enforcement of an inclined boundary condition
- ADPressureApplies a pressure on a given boundary in a given direction
- ADTorqueApply a moment as tractions distributed over a surface around a pivot point. This should operate on the displaced mesh for large deformations.
- CoupledPressureBCApplies a pressure from a variable on a given boundary in a given direction
- DashpotBCModel a dashpot boundary condition where the traction is proportional to the normal velocity.
- DirectDirichletBCImposes the essential boundary condition , where is a constant, controllable value.
- DirectFunctionDirichletBCImposes the essential boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- DisplacementAboutAxisImplements a boundary condition that enforces rotationaldisplacement around an axis on a boundary
- ExplicitDirichletBCImposes the essential boundary condition , where is a constant, controllable value.
- ExplicitFunctionDirichletBCImposes the essential boundary condition , where is a (possibly) time and space-dependent MOOSE Function.
- InteractionIntegralBenchmarkBCImplements a boundary condition that enforces a displacement field around a crack tip based on applied stress intensity factors.
- PenaltyInclinedNoDisplacementBCPenalty Enforcement of an inclined boundary condition
- PresetAccelerationPrescribe acceleration on a given boundary in a given direction
- PresetDisplacementPrescribe the displacement on a given boundary in a given direction.
- PresetVelocitySets the boundary displacements through time from an imposed velocity
- PressureApplies a pressure on a given boundary in a given direction
- StickyBCImposes the boundary condition if exceeds the bounds provided
- TorqueApply a moment as tractions distributed over a surface around a pivot point. This should operate on the displaced mesh for large deformations.
- Thermal Hydraulics App
- ADBoundaryFlux3EqnBCBoundary conditions for the 1-D, 1-phase, variable-area Euler equations
- ADConvectionHeatTransfer3DBCAdds a convective heat flux boundary condition between the local component heat structure and a 3D heat structure
- ADConvectionHeatTransferBCAdds a convective heat flux boundary condition with user-specified ambient temperature and heat transfer coefficient functions
- ADConvectionHeatTransferRZBCConvection BC for RZ domain in XY coordinate system
- ADExternalAppConvectionHeatTransferBCConvection BC from an external application
- ADExternalAppConvectionHeatTransferRZBCConvection BC from an external application for RZ domain in XY coordinate system
- ADGateValve1PhaseBCAdds boundary fluxes for flow channels connected to a 1-phase gate valve
- ADHSHeatFluxBCApplies a specified heat flux to the side of a plate heat structure
- ADHSHeatFluxRZBCApplies a specified heat flux to the side of a cylindrical heat structure in XY coordinates
- ADHeatFlux3EqnBCWall heat flux boundary condition for the energy equation
- ADJunctionOneToOne1PhaseBCAdds boundary fluxes for flow channels connected to a 1-phase one-to-one junction
- ADRadiativeHeatFluxBCRadiative heat transfer boundary condition for a plate heat structure
- ADRadiativeHeatFluxRZBCRadiative heat transfer boundary condition for a cylindrical heat structure
- ADVolumeJunction1PhaseBCAdds boundary fluxes for flow channels connected to a 1-phase volume junction
- BoundaryFlux3EqnBCBoundary conditions for the 1-D, 1-phase, variable-area Euler equations
- BoundaryFluxGasMixBCBoundary conditions for a FlowChannelGasMix using a boundary flux object.
- ConvectionHeatTransferBCAdds a convective heat flux boundary condition with user-specified ambient temperature and heat transfer coefficient functions
- ConvectionHeatTransferRZBCConvection BC for RZ domain in XY coordinate system
- ExternalAppConvectionHeatTransferBCConvection BC from an external application
- ExternalAppConvectionHeatTransferRZBCConvection BC from an external application for RZ domain in XY coordinate system
- HSCoupler2D2DRadiationRZBCAdds boundary heat flux terms for HSCoupler2D2DRadiation
- HSCoupler2D3DBCAdds boundary heat flux terms for HSCoupler2D3D
- HeatStructure2DCouplerBCApplies BC for HeatStructure2DCoupler for plate heat structure
- HeatStructure2DCouplerRZBCApplies BC for HeatStructure2DCoupler for cylindrical heat structure in a XY coordinate system
- HeatStructure2DRadiationCouplerRZBCApplies BC for HeatStructure2DRadiationCouplerRZ
- KokkosRadiativeHeatFluxBCRadiative heat transfer boundary condition for a plate heat structure
- RadiativeHeatFluxBCRadiative heat transfer boundary condition for a plate heat structure
- RadiativeHeatFluxRZBCRadiative heat transfer boundary condition for a cylindrical heat structure in a XY coordinate system
- Chemical Reactions App
- ChemicalOutFlowBCChemical flux boundary condition
- Phase Field App
- ADPhaseFieldContactAngleBCEnforce contact angle BC using phase field variable