Implements all the methods for assembling a hybridized local discontinuous Galerkin (LDG-H), which is a type of HDG method, discretization of the diffusion equation. More...
#include <DiffusionLHDGAssemblyHelper.h>
Public Member Functions | |
DiffusionLHDGAssemblyHelper (const MooseObject *const moose_obj, MaterialPropertyInterface *const mpi, MooseVariableDependencyInterface *const mvdi, const TransientInterface *const ti, const FEProblemBase &fe_problem, SystemBase &sys, const THREAD_ID tid) | |
void | checkCoupling () |
Static Public Member Functions | |
static InputParameters | validParams () |
static std::string | deduceFunctorName (const std::string &name, const InputParameters ¶ms) |
Helper to look up a functor name through the input parameter keys. More... | |
Protected Member Functions | |
void | vectorVolumeResidual (const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseVector< Number > &vector_re) |
Computes a local residual vector for the weak form: (q, v) + (u, div(v)) where q is the vector field representing the gradient of u, v are its associated test functions, and u is the diffused scalar field. More... | |
void | vectorVolumeJacobian (const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &vector_vector_jac, DenseMatrix< Number > &vector_scalar_jac) |
Computes a local Jacobian matrix for the weak form: (q, v) + (u, div(v)) where q is the vector field representing the gradient, v are its associated test functions, and u is the scalar field. More... | |
void | scalarVolumeResidual (const MooseArray< Gradient > &vector_field, const Moose::Functor< Real > &source, const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const Elem *const current_elem, const MooseArray< Point > &q_point, DenseVector< Number > &scalar_re) |
Computes a local residual vector for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity, w are the test functions associated with the scalar field, and f is a forcing function. More... | |
void | scalarVolumeJacobian (const MooseArray< Real > &JxW, const libMesh::QBase &qrule, DenseMatrix< Number > &scalar_vector_jac) |
Computes a local Jacobian matrix for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity, w are the test functions associated with the scalar field, and f is a forcing function. More... | |
void | vectorFaceResidual (const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &vector_re) |
Computes a local residual vector for the weak form: -<{u}, n*v> where {u} is the trace of the scalar field, n is the normal vector, and v are the test functions associated with the gradient field. More... | |
void | vectorFaceJacobian (const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &vector_lm_jac) |
Computes a local Jacobian matrix for the weak form: -<{u}, n*v> where {u} is the trace of the scalar field, n is the normal vector, and v are the test functions associated with the gradient field. More... | |
void | scalarFaceResidual (const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &scalar_re) |
Computes a local residual vector for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n, w> More... | |
void | scalarFaceJacobian (const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &scalar_vector_jac, DenseMatrix< Number > &scalar_scalar_jac, DenseMatrix< Number > &scalar_lm_jac) |
Computes a local Jacobian matrix for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n, w> More... | |
void | lmFaceResidual (const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const MooseArray< Number > &lm_sol, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseVector< Number > &lm_re) |
Computes a local residual vector for the weak form: -<Dq*n, > + < * (u - {u}) * n * n, > More... | |
void | lmFaceJacobian (const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &lm_vec_jac, DenseMatrix< Number > &lm_scalar_jac, DenseMatrix< Number > &lm_lm_jac) |
Computes a local Jacobian matrix for the weak form: -<Dq*n, > + < * (u - {u}) * n * n, > More... | |
void | vectorDirichletResidual (const Moose::Functor< Real > &dirichlet_value, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const current_elem, const unsigned int current_side, const MooseArray< Point > &q_point_face, DenseVector< Number > &vector_re) |
Weakly imposes a Dirichlet condition for the scalar field in the vector (gradient) equation. More... | |
void | scalarDirichletResidual (const MooseArray< Gradient > &vector_sol, const MooseArray< Number > &scalar_sol, const Moose::Functor< Real > &dirichlet_value, const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, const Elem *const current_elem, const unsigned int current_side, const MooseArray< Point > &q_point_face, DenseVector< Number > &scalar_re) |
Weakly imposes a Dirichlet condition for the scalar field in the scalar field equation. More... | |
void | scalarDirichletJacobian (const MooseArray< Real > &JxW_face, const libMesh::QBase &qrule_face, const MooseArray< Point > &normals, DenseMatrix< Number > &scalar_vector_jac, DenseMatrix< Number > &scalar_scalar_jac) |
Computes the Jacobian for a Dirichlet condition for the scalar field in the scalar field equation. More... | |
void | createIdentityResidual (const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const MooseArray< std::vector< Real >> &phi, const MooseArray< Number > &sol, DenseVector< Number > &re) |
Creates residuals corresponding to the weak form (v, {u}), or stated simply this routine can be used to drive Lagrange multiplier values on the boundary to zero. More... | |
void | createIdentityJacobian (const MooseArray< Real > &JxW, const libMesh::QBase &qrule, const MooseArray< std::vector< Real >> &phi, DenseMatrix< Number > &ke) |
As above, but for the Jacobians. More... | |
std::string | deduceFunctorName (const std::string &name) const |
Small helper to look up a functor name through the input parameter keys. More... | |
template<typename T > | |
const Moose::Functor< T > & | getFunctor (const std::string &name) |
Retrieves a functor from the subproblem. More... | |
template<typename T > | |
const Moose::Functor< T > & | getFunctor (const std::string &name, THREAD_ID tid) |
Retrieves a functor from the subproblem. More... | |
template<typename T > | |
const Moose::Functor< T > & | getFunctor (const std::string &name, SubProblem &subproblem) |
Retrieves a functor from the passed-in subproblem. More... | |
template<typename T > | |
const Moose::Functor< T > & | getFunctor (const std::string &name, SubProblem &subproblem, THREAD_ID tid) |
Retrieves a functor from the passed-in subproblem. More... | |
bool | isFunctor (const std::string &name) const |
Checks the subproblem for the given functor. More... | |
bool | isFunctor (const std::string &name, const SubProblem &subproblem) const |
Checks the passed-in subproblem for the given functor. More... | |
Moose::ElemArg | makeElemArg (const Elem *elem, bool correct_skewnewss=false) const |
Helper method to create an elemental argument for a functor that includes whether to perform skewness corrections. More... | |
template<typename T > | |
void | checkFunctorSupportsSideIntegration (const std::string &name, bool qp_integration) |
Throws error if the functor does not support the requested side integration. More... | |
Private Attributes | |
const MooseObject & | _moose_obj |
A reference to our associated MooseObject for error reporting. More... | |
const FEProblemBase & | _dhah_fe_problem |
A reference to the finite element problem used for coupling checks. More... | |
const SystemBase & | _dhah_sys |
A reference to the nonlinear system used for coupling checks. More... | |
Implements all the methods for assembling a hybridized local discontinuous Galerkin (LDG-H), which is a type of HDG method, discretization of the diffusion equation.
These routines may be called by both HDG kernels and integrated boundary conditions. The implementation here is based (but not exactly based) on "A superconvergent LDG-hybridizable Galerkin method for second-order elliptic problems" by Cockburn
Definition at line 40 of file DiffusionLHDGAssemblyHelper.h.
DiffusionLHDGAssemblyHelper::DiffusionLHDGAssemblyHelper | ( | const MooseObject *const | moose_obj, |
MaterialPropertyInterface *const | mpi, | ||
MooseVariableDependencyInterface *const | mvdi, | ||
const TransientInterface *const | ti, | ||
const FEProblemBase & | fe_problem, | ||
SystemBase & | sys, | ||
const THREAD_ID | tid | ||
) |
Definition at line 37 of file DiffusionLHDGAssemblyHelper.C.
void DiffusionLHDGAssemblyHelper::checkCoupling | ( | ) |
Definition at line 77 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::initialSetup(), DiffusionLHDGPrescribedGradientBC::initialSetup(), and DiffusionLHDGKernel::initialSetup().
|
protectedinherited |
Throws error if the functor does not support the requested side integration.
[in] | name | Name of functor or functor parameter |
[in] | qp_integration | True if performing qp integration, false if face info |
Definition at line 236 of file FunctorInterface.h.
|
protected |
As above, but for the Jacobians.
Definition at line 418 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::computeJacobian().
|
protected |
Creates residuals corresponding to the weak form (v, {u}), or stated simply this routine can be used to drive Lagrange multiplier values on the boundary to zero.
This should be used on boundaries where there are Dirichlet conditions for the primal variables such that there is no need for the Lagrange multiplier variables
Definition at line 403 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::computeResidual().
|
staticinherited |
Helper to look up a functor name through the input parameter keys.
name | The input parameter name that we are trying to deduce the functor name for |
params | The input parameters object that we will be checking for parameters named name |
Definition at line 28 of file FunctorInterface.C.
Referenced by FunctorInterface::checkFunctorSupportsSideIntegration(), FunctorInterface::deduceFunctorName(), FunctorInterface::getFunctor(), and FunctorInterface::isFunctor().
|
protectedinherited |
Small helper to look up a functor name through the input parameter keys.
Definition at line 60 of file FunctorInterface.C.
|
protectedinherited |
Retrieves a functor from the subproblem.
This method also leverages the ability to create default functors if the user passed an integer or real in the input file
name | The name of the functor to retrieve. This should match the functor parameter name, not the actual name of the functor created in the input file |
Definition at line 200 of file FunctorInterface.h.
Referenced by MaterialFunctorConverterTempl< T >::MaterialFunctorConverterTempl().
|
protectedinherited |
Retrieves a functor from the subproblem.
This method also leverages the ability to create default functors if the user passed an integer or real in the input file
name | The name of the functor to retrieve. This should match the functor parameter name, not the actual name of the functor created in the input file |
tid | The thread ID used to retrieve the functor from this interface's subproblem |
Definition at line 192 of file FunctorInterface.h.
|
protectedinherited |
Retrieves a functor from the passed-in subproblem.
This method also leverages the ability to create default functors if the user passed an integer or real in the input file
name | The name of the functor to retrieve. This should match the functor parameter name, not the actual name of the functor created in the input file |
subproblem | The subproblem to query for the functor |
Definition at line 185 of file FunctorInterface.h.
|
protectedinherited |
Retrieves a functor from the passed-in subproblem.
This method also leverages the ability to create default functors if the user passed an integer or real in the input file
name | The name of the functor to retrieve. This should match the functor parameter name, not the actual name of the functor created in the input file |
subproblem | The subproblem to query for the functor |
tid | The thread ID used to retrieve the functor from the subproblem |
Definition at line 176 of file FunctorInterface.h.
|
protectedinherited |
Checks the subproblem for the given functor.
This will not query default functors potentially stored in this object, e.g. this method will return false if the user passed an int or real to the functor param in the input file
name | The name of the functor to check. This should match the functor parameter name, not the actual name of the functor created in the input file |
Definition at line 113 of file FunctorInterface.C.
|
protectedinherited |
Checks the passed-in subproblem for the given functor.
This will not query default functors potentially stored in this object, e.g. this method will return false if the user passed an int or real to the functor param in the input file
name | The name of the functor to check. This should match the functor parameter name, not the actual name of the functor created in the input file |
subproblem | The subproblem to query for the functor |
Definition at line 104 of file FunctorInterface.C.
|
protected |
Computes a local Jacobian matrix for the weak form: -<Dq*n, > + < * (u - {u}) * n * n, >
Definition at line 301 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Computes a local residual vector for the weak form: -<Dq*n, > + < * (u - {u}) * n * n, >
Definition at line 277 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protectedinherited |
Helper method to create an elemental argument for a functor that includes whether to perform skewness corrections.
Definition at line 120 of file FunctorInterface.C.
Referenced by LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryNormalGradient(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValue(), LinearFVReaction::computeMatrixContribution(), LinearFVTimeDerivative::computeMatrixContribution(), FVFunctorTimeKernel::computeQpResidual(), FVCoupledForce::computeQpResidual(), FVMassMatrix::computeQpResidual(), FVIntegralValueConstraint::computeQpResidual(), FVBoundedValueConstraint::computeQpResidual(), FVPointValueConstraint::computeQpResidual(), LinearFVSource::computeRightHandSideContribution(), SecondTimeDerivativeAux::computeValue(), TimeDerivativeAux::computeValue(), FunctorAux::computeValue(), FunctorCoordinatesFunctionAux::computeValue(), FunctorTimes::initialize(), and LinearFVTimeDerivative::setCurrentElemInfo().
|
protected |
Computes the Jacobian for a Dirichlet condition for the scalar field in the scalar field equation.
Definition at line 379 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::computeJacobian().
|
protected |
Weakly imposes a Dirichlet condition for the scalar field in the scalar field equation.
Definition at line 352 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::computeResidual().
|
protected |
Computes a local Jacobian matrix for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n, w>
Definition at line 249 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Computes a local residual vector for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n, w>
Definition at line 225 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Computes a local Jacobian matrix for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity, w are the test functions associated with the scalar field, and f is a forcing function.
Definition at line 172 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGKernel::computeJacobian().
|
protected |
Computes a local residual vector for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity, w are the test functions associated with the scalar field, and f is a forcing function.
Definition at line 145 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGKernel::computeResidual().
|
static |
Definition at line 22 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::validParams(), DiffusionLHDGKernel::validParams(), and DiffusionLHDGPrescribedGradientBC::validParams().
|
protected |
Weakly imposes a Dirichlet condition for the scalar field in the vector (gradient) equation.
Definition at line 329 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGDirichletBC::computeResidual().
|
protected |
Computes a local Jacobian matrix for the weak form: -<{u}, n*v> where {u} is the trace of the scalar field, n is the normal vector, and v are the test functions associated with the gradient field.
Definition at line 206 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Computes a local residual vector for the weak form: -<{u}, n*v> where {u} is the trace of the scalar field, n is the normal vector, and v are the test functions associated with the gradient field.
Definition at line 190 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Computes a local Jacobian matrix for the weak form: (q, v) + (u, div(v)) where q is the vector field representing the gradient, v are its associated test functions, and u is the scalar field.
Definition at line 124 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGKernel::computeJacobian().
|
protected |
Computes a local residual vector for the weak form: (q, v) + (u, div(v)) where q is the vector field representing the gradient of u, v are its associated test functions, and u is the diffused scalar field.
Definition at line 102 of file DiffusionLHDGAssemblyHelper.C.
Referenced by DiffusionLHDGKernel::computeResidual().
|
protected |
A data member used for determining when to compute the Jacobian.
Definition at line 266 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeOffDiagJacobian(), DiffusionLHDGPrescribedGradientBC::computeOffDiagJacobian(), DiffusionLHDGKernel::computeOffDiagJacobian(), DiffusionLHDGDirichletBC::jacobianSetup(), DiffusionLHDGPrescribedGradientBC::jacobianSetup(), and DiffusionLHDGKernel::jacobianSetup().
|
private |
A reference to the finite element problem used for coupling checks.
Definition at line 281 of file DiffusionLHDGAssemblyHelper.h.
Referenced by checkCoupling().
|
private |
A reference to the nonlinear system used for coupling checks.
Definition at line 284 of file DiffusionLHDGAssemblyHelper.h.
Referenced by checkCoupling().
|
protected |
The diffusivity.
Definition at line 257 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeResidual(), lmFaceJacobian(), lmFaceResidual(), scalarDirichletJacobian(), scalarDirichletResidual(), scalarFaceJacobian(), scalarFaceResidual(), scalarVolumeJacobian(), and scalarVolumeResidual().
|
protected |
Definition at line 249 of file DiffusionLHDGAssemblyHelper.h.
Referenced by vectorVolumeJacobian(), and vectorVolumeResidual().
|
protected |
Definition at line 248 of file DiffusionLHDGAssemblyHelper.h.
Referenced by scalarVolumeJacobian(), and scalarVolumeResidual().
|
protected |
Definition at line 232 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::additionalROVariables(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), DiffusionLHDGKernel::computeResidualOnSide(), and DiffusionLHDGAssemblyHelper().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Definition at line 254 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), lmFaceJacobian(), lmFaceResidual(), scalarFaceJacobian(), and vectorFaceJacobian().
|
protected |
Definition at line 269 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Definition at line 238 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 243 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
private |
A reference to our associated MooseObject for error reporting.
Definition at line 278 of file DiffusionLHDGAssemblyHelper.h.
Referenced by checkCoupling().
|
protected |
Definition at line 236 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), DiffusionLHDGKernel::computeResidualOnSide(), scalarDirichletJacobian(), and vectorDirichletResidual().
|
protected |
Definition at line 241 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Definition at line 247 of file DiffusionLHDGAssemblyHelper.h.
Referenced by scalarVolumeResidual(), and vectorVolumeJacobian().
|
protected |
Definition at line 253 of file DiffusionLHDGAssemblyHelper.h.
Referenced by lmFaceJacobian(), scalarDirichletJacobian(), scalarDirichletResidual(), scalarFaceJacobian(), and scalarFaceResidual().
|
protected |
Definition at line 269 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Our stabilization coefficient.
Definition at line 263 of file DiffusionLHDGAssemblyHelper.h.
Referenced by lmFaceJacobian(), lmFaceResidual(), scalarDirichletJacobian(), scalarDirichletResidual(), scalarFaceJacobian(), and scalarFaceResidual().
|
protected |
Reference to transient interface.
Definition at line 260 of file DiffusionLHDGAssemblyHelper.h.
Referenced by scalarDirichletResidual(), scalarVolumeResidual(), and vectorDirichletResidual().
|
protected |
Definition at line 237 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), DiffusionLHDGKernel::computeResidualOnSide(), scalarDirichletJacobian(), and scalarDirichletResidual().
|
protected |
Definition at line 233 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::additionalROVariables(), DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), DiffusionLHDGKernel::computeResidualOnSide(), and DiffusionLHDGAssemblyHelper().
|
protected |
Definition at line 242 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 231 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGDirichletBC::computeJacobian(), DiffusionLHDGPrescribedGradientBC::computeJacobian(), DiffusionLHDGKernel::computeJacobian(), DiffusionLHDGKernel::computeJacobianOnSide(), DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobianOnSide().
|
protected |
Definition at line 246 of file DiffusionLHDGAssemblyHelper.h.
Referenced by scalarVolumeJacobian(), vectorVolumeJacobian(), and vectorVolumeResidual().
|
protected |
Definition at line 252 of file DiffusionLHDGAssemblyHelper.h.
Referenced by lmFaceJacobian(), scalarDirichletJacobian(), scalarFaceJacobian(), vectorDirichletResidual(), vectorFaceJacobian(), and vectorFaceResidual().
|
protected |
Definition at line 269 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGKernel::computeResidual(), DiffusionLHDGDirichletBC::computeResidual(), DiffusionLHDGPrescribedGradientBC::computeResidual(), and DiffusionLHDGKernel::computeResidualOnSide().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobian().
|
protected |
Definition at line 272 of file DiffusionLHDGAssemblyHelper.h.
Referenced by DiffusionLHDGPrescribedGradientBC::computeJacobian(), and DiffusionLHDGKernel::computeJacobian().