https://mooseframework.inl.gov
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
DiffusionLHDGAssemblyHelper Class Reference

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>

Inheritance diagram for DiffusionLHDGAssemblyHelper:
[legend]

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 &params)
 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...
 

Protected Attributes

const MooseVariableFE< Real > & _u_var
 
const MooseVariableFE< RealVectorValue > & _grad_u_var
 
const MooseVariableFE< Real > & _u_face_var
 
const std::vector< dof_id_type > & _qu_dof_indices
 
const std::vector< dof_id_type > & _u_dof_indices
 
const std::vector< dof_id_type > & _lm_u_dof_indices
 
const MooseArray< libMesh::Gradient > & _qu_sol
 
const MooseArray< Number > & _u_sol
 
const MooseArray< Number > & _lm_u_sol
 
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
 
const MooseArray< std::vector< Real > > & _scalar_phi
 
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
 
const MooseArray< std::vector< Real > > & _div_vector_phi
 
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
 
const MooseArray< std::vector< Real > > & _scalar_phi_face
 
const MooseArray< std::vector< Real > > & _lm_phi_face
 
const MaterialProperty< Real > & _diff
 The diffusivity. More...
 
const TransientInterface_ti
 Reference to transient interface. More...
 
const Real _tau
 Our stabilization coefficient. More...
 
const Elem * _cached_elem
 A data member used for determining when to compute the Jacobian. More...
 
DenseVector< Number_vector_re
 
DenseVector< Number_scalar_re
 
DenseVector< Number_lm_re
 
DenseMatrix< Number_vector_vector_jac
 
DenseMatrix< Number_vector_scalar_jac
 
DenseMatrix< Number_scalar_vector_jac
 
DenseMatrix< Number_scalar_scalar_jac
 
DenseMatrix< Number_scalar_lm_jac
 
DenseMatrix< Number_lm_scalar_jac
 
DenseMatrix< Number_lm_lm_jac
 
DenseMatrix< Number_vector_lm_jac
 
DenseMatrix< Number_lm_vector_jac
 

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...
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ DiffusionLHDGAssemblyHelper()

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.

45  : ADFunctorInterface(moose_obj),
46  _u_var(sys.getFieldVariable<Real>(tid, moose_obj->getParam<NonlinearVariableName>("variable"))),
48  tid, moose_obj->getParam<NonlinearVariableName>("gradient_variable"))),
49  _u_face_var(sys.getFieldVariable<Real>(
50  tid, moose_obj->getParam<NonlinearVariableName>("face_variable"))),
55  _u_sol(_u_var.sln()),
64  _diff(mpi->getMaterialProperty<Real>("diffusivity")),
65  _ti(*ti),
66  _tau(moose_obj->getParam<Real>("tau")),
67  _cached_elem(nullptr),
68  _moose_obj(*moose_obj),
69  _dhah_fe_problem(fe_problem),
70  _dhah_sys(sys)
71 {
74 }
const FEProblemBase & _dhah_fe_problem
A reference to the finite element problem used for coupling checks.
const std::vector< dof_id_type > & _lm_u_dof_indices
const MooseArray< std::vector< Real > > & _scalar_phi_face
const MooseArray< Number > & _u_sol
const MooseArray< std::vector< Real > > & _scalar_phi
const MooseArray< std::vector< Real > > & _div_vector_phi
const std::vector< dof_id_type > & _u_dof_indices
const Elem * _cached_elem
A data member used for determining when to compute the Jacobian.
const MooseArray< std::vector< Real > > & _lm_phi_face
const MooseArray< libMesh::Gradient > & _qu_sol
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
const FieldVariablePhiDivergence & divPhi() const override final
Divergence of the shape functions.
const std::vector< dof_id_type > & _qu_dof_indices
const SystemBase & _dhah_sys
A reference to the nonlinear system used for coupling checks.
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
const FieldVariableValue & sln() const override
element solutions
const FieldVariablePhiValue & phi() const override
Return the variable&#39;s elemental shape functions.
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
const MaterialProperty< Real > & _diff
The diffusivity.
const TransientInterface & _ti
Reference to transient interface.
const MooseVariableFE< Real > & _u_face_var
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
const Real _tau
Our stabilization coefficient.
MooseVariableFE< T > & getFieldVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:110
const MooseVariableFE< Real > & _u_var
const MooseObject & _moose_obj
A reference to our associated MooseObject for error reporting.
const MooseArray< Number > & _lm_u_sol
ADFunctorInterface(const MooseObject *moose_object)
const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
const MooseVariableFE< RealVectorValue > & _grad_u_var
const FieldVariablePhiGradient & gradPhi() const override final
Return the gradients of the variable&#39;s elemental shape functions.
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
const MaterialProperty< T > & getMaterialProperty(const std::string &name, const unsigned int state=0)

Member Function Documentation

◆ checkCoupling()

void DiffusionLHDGAssemblyHelper::checkCoupling ( )

Definition at line 77 of file DiffusionLHDGAssemblyHelper.C.

Referenced by DiffusionLHDGDirichletBC::initialSetup(), DiffusionLHDGPrescribedGradientBC::initialSetup(), and DiffusionLHDGKernel::initialSetup().

78 {
79  // This check must occur after FEProblemBase::init()
81  return;
83  {
84  const auto * const cm = _dhah_fe_problem.couplingMatrix(_dhah_sys.number());
85  for (const auto i : make_range(cm->size()))
86  for (const auto j : make_range(cm->size()))
87  if ((*cm)(i, j) != true)
88  goto error;
89 
90  return;
91  }
92 
93 error:
95  "This class encodes the full Jacobian regardless of user input file specification, "
96  "so please request full coupling for system ",
97  _dhah_sys.name(),
98  " in your Preconditioning block for consistency");
99 }
const FEProblemBase & _dhah_fe_problem
A reference to the finite element problem used for coupling checks.
const libMesh::CouplingMatrix * couplingMatrix(const unsigned int nl_sys_num) const override
The coupling matrix defining what blocks exist in the preconditioning matrix.
virtual const std::string & name() const
Definition: SystemBase.C:1330
const SystemBase & _dhah_sys
A reference to the nonlinear system used for coupling checks.
Moose::CouplingType coupling() const
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1159
const MooseObject & _moose_obj
A reference to our associated MooseObject for error reporting.
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.

◆ checkFunctorSupportsSideIntegration()

template<typename T >
void FunctorInterface::checkFunctorSupportsSideIntegration ( const std::string &  name,
bool  qp_integration 
)
protectedinherited

Throws error if the functor does not support the requested side integration.

Parameters
[in]nameName of functor or functor parameter
[in]qp_integrationTrue if performing qp integration, false if face info

Definition at line 236 of file FunctorInterface.h.

237 {
238  const std::string functor_name = deduceFunctorName(name);
239  const auto & functor = getFunctor<T>(name);
240  if (qp_integration)
241  {
242  if (!functor.supportsElemSideQpArg())
243  mooseError("Quadrature point integration was requested, but the functor '",
244  functor_name,
245  "' does not support this.");
246  }
247  else
248  {
249  if (!functor.supportsFaceArg())
250  mooseError("Face info integration was requested, but the functor '",
251  functor_name,
252  "' does not support this.");
253  }
254 }
std::string name(const ElemQuality q)
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
static std::string deduceFunctorName(const std::string &name, const InputParameters &params)
Helper to look up a functor name through the input parameter keys.

◆ createIdentityJacobian()

void DiffusionLHDGAssemblyHelper::createIdentityJacobian ( const MooseArray< Real > &  JxW,
const libMesh::QBase qrule,
const MooseArray< std::vector< Real >> &  phi,
DenseMatrix< Number > &  ke 
)
protected

As above, but for the Jacobians.

Definition at line 418 of file DiffusionLHDGAssemblyHelper.C.

Referenced by DiffusionLHDGDirichletBC::computeJacobian().

422 {
423  for (const auto qp : make_range(qrule.n_points()))
424  for (const auto i : index_range(phi))
425  {
426  const auto qpi_term = JxW[qp] * phi[i][qp];
427  for (const auto j : index_range(phi))
428  ke(i, j) -= phi[j][qp] * qpi_term;
429  }
430 }
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)

◆ createIdentityResidual()

void DiffusionLHDGAssemblyHelper::createIdentityResidual ( const MooseArray< Real > &  JxW,
const libMesh::QBase qrule,
const MooseArray< std::vector< Real >> &  phi,
const MooseArray< Number > &  sol,
DenseVector< Number > &  re 
)
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().

408 {
409  for (const auto qp : make_range(qrule.n_points()))
410  {
411  const auto qp_term = JxW[qp] * sol[qp];
412  for (const auto i : index_range(phi))
413  re(i) -= phi[i][qp] * qp_term;
414  }
415 }
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)

◆ deduceFunctorName() [1/2]

std::string FunctorInterface::deduceFunctorName ( const std::string &  name,
const InputParameters params 
)
staticinherited

Helper to look up a functor name through the input parameter keys.

Parameters
nameThe input parameter name that we are trying to deduce the functor name for
paramsThe input parameters object that we will be checking for parameters named name
Returns
The functor name

Definition at line 28 of file FunctorInterface.C.

Referenced by FunctorInterface::checkFunctorSupportsSideIntegration(), FunctorInterface::deduceFunctorName(), FunctorInterface::getFunctor(), and FunctorInterface::isFunctor().

29 {
30  if (params.isParamValid(name))
31  {
32  if (params.have_parameter<MooseFunctorName>(name))
33  return params.get<MooseFunctorName>(name);
34  // variables, functor material properties, functions, and post-processors are also functors
35  else if (params.have_parameter<MaterialPropertyName>(name))
36  return params.get<MaterialPropertyName>(name);
37  else if (params.have_parameter<VariableName>(name))
38  return params.get<VariableName>(name);
39  else if (params.have_parameter<std::vector<VariableName>>(name))
40  {
41  const auto & var_names = params.get<std::vector<VariableName>>(name);
42  if (var_names.size() != 1)
43  mooseError("We only support a single variable name for retrieving a functor");
44  return var_names[0];
45  }
46  else if (params.have_parameter<NonlinearVariableName>(name))
47  return params.get<NonlinearVariableName>(name);
48  else if (params.have_parameter<FunctionName>(name))
49  return params.get<FunctionName>(name);
50  else if (params.have_parameter<PostprocessorName>(name))
51  return params.get<PostprocessorName>(name);
52  else
53  mooseError("Invalid parameter type for retrieving a functor");
54  }
55  else
56  return name;
57 }
std::string name(const ElemQuality q)
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
bool have_parameter(std::string_view name) const
A wrapper around the Parameters base class method.
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.

◆ deduceFunctorName() [2/2]

std::string FunctorInterface::deduceFunctorName ( const std::string &  name) const
protectedinherited

Small helper to look up a functor name through the input parameter keys.

Definition at line 60 of file FunctorInterface.C.

61 {
62  return deduceFunctorName(name, _fi_params);
63 }
const InputParameters & _fi_params
Parameters of the object with this interface.
static std::string deduceFunctorName(const std::string &name, const InputParameters &params)
Helper to look up a functor name through the input parameter keys.

◆ getFunctor() [1/4]

template<typename T >
const Moose::Functor< T > & FunctorInterface::getFunctor ( const std::string &  name)
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

Parameters
nameThe 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
Returns
The functor

Definition at line 200 of file FunctorInterface.h.

Referenced by MaterialFunctorConverterTempl< T >::MaterialFunctorConverterTempl().

201 {
202  mooseAssert(_fi_subproblem, "This must be non-null");
203  return getFunctor<T>(name, *_fi_subproblem, _fi_tid);
204 }
std::string name(const ElemQuality q)
SubProblem *const _fi_subproblem
Pointer to subproblem if the subproblem pointer parameter was set.
const THREAD_ID _fi_tid
Current threaded it.

◆ getFunctor() [2/4]

template<typename T >
const Moose::Functor< T > & FunctorInterface::getFunctor ( const std::string &  name,
THREAD_ID  tid 
)
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

Parameters
nameThe 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
tidThe thread ID used to retrieve the functor from this interface's subproblem
Returns
The functor

Definition at line 192 of file FunctorInterface.h.

193 {
194  mooseAssert(_fi_subproblem, "This must be non-null");
195  return getFunctor<T>(name, *_fi_subproblem, tid);
196 }
std::string name(const ElemQuality q)
SubProblem *const _fi_subproblem
Pointer to subproblem if the subproblem pointer parameter was set.

◆ getFunctor() [3/4]

template<typename T >
const Moose::Functor< T > & FunctorInterface::getFunctor ( const std::string &  name,
SubProblem subproblem 
)
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

Parameters
nameThe 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
subproblemThe subproblem to query for the functor
Returns
The functor

Definition at line 185 of file FunctorInterface.h.

186 {
187  return getFunctor<T>(name, subproblem, _fi_tid);
188 }
std::string name(const ElemQuality q)
const THREAD_ID _fi_tid
Current threaded it.

◆ getFunctor() [4/4]

template<typename T >
const Moose::Functor< T > & FunctorInterface::getFunctor ( const std::string &  name,
SubProblem subproblem,
THREAD_ID  tid 
)
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

Parameters
nameThe 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
subproblemThe subproblem to query for the functor
tidThe thread ID used to retrieve the functor from the subproblem
Returns
The functor

Definition at line 176 of file FunctorInterface.h.

177 {
178  // Check if the supplied parameter is a valid input parameter key
179  std::string functor_name = deduceFunctorName(name);
180  return getFunctorByName<T>(functor_name, subproblem, tid);
181 }
static std::string deduceFunctorName(const std::string &name, const InputParameters &params)
Helper to look up a functor name through the input parameter keys.

◆ isFunctor() [1/2]

bool FunctorInterface::isFunctor ( const std::string &  name) const
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

Parameters
nameThe 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
Returns
Whether the subproblem has the specified functor

Definition at line 113 of file FunctorInterface.C.

114 {
115  mooseAssert(_fi_subproblem, "This must be non-null");
116  return isFunctor(name, *_fi_subproblem);
117 }
SubProblem *const _fi_subproblem
Pointer to subproblem if the subproblem pointer parameter was set.
bool isFunctor(const std::string &name) const
Checks the subproblem for the given functor.

◆ isFunctor() [2/2]

bool FunctorInterface::isFunctor ( const std::string &  name,
const SubProblem subproblem 
) const
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

Parameters
nameThe 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
subproblemThe subproblem to query for the functor
Returns
Whether the subproblem has the specified functor

Definition at line 104 of file FunctorInterface.C.

105 {
106  // Check if the supplied parameter is a valid input parameter key
107  std::string functor_name = deduceFunctorName(name);
108 
109  return subproblem.hasFunctor(functor_name, _fi_tid);
110 }
bool hasFunctor(const std::string &name, const THREAD_ID tid) const
checks whether we have a functor corresponding to name on the thread id tid
Definition: SubProblem.C:1264
const THREAD_ID _fi_tid
Current threaded it.
static std::string deduceFunctorName(const std::string &name, const InputParameters &params)
Helper to look up a functor name through the input parameter keys.

◆ lmFaceJacobian()

void DiffusionLHDGAssemblyHelper::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 
)
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().

307 {
308  for (const auto qp : make_range(qrule_face.n_points()))
309  {
310  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * normals[qp];
311  const auto stab_qp_term = JxW_face[qp] * _tau * normals[qp] * normals[qp];
312 
313  for (const auto i : make_range(lm_vec_jac.m()))
314  {
315  const auto vector_qpi_term = vector_qp_term * _lm_phi_face[i][qp];
316  for (const auto j : make_range(lm_vec_jac.n()))
317  lm_vec_jac(i, j) -= vector_qpi_term * _vector_phi_face[j][qp];
318 
319  const auto lm_qpi_term = stab_qp_term * _lm_phi_face[i][qp];
320  for (const auto j : make_range(lm_scalar_jac.n()))
321  lm_scalar_jac(i, j) += lm_qpi_term * _scalar_phi_face[j][qp];
322  for (const auto j : make_range(lm_lm_jac.n()))
323  lm_lm_jac(i, j) -= lm_qpi_term * _lm_phi_face[j][qp];
324  }
325  }
326 }
const MooseArray< std::vector< Real > > & _scalar_phi_face
const MooseArray< std::vector< Real > > & _lm_phi_face
unsigned int m() const
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
unsigned int n() const
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face

◆ lmFaceResidual()

void DiffusionLHDGAssemblyHelper::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 
)
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().

284 {
285  for (const auto qp : make_range(qrule_face.n_points()))
286  {
287  // vector
288  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * (vector_sol[qp] * normals[qp]);
289  // stabilization term
290  const auto stab_qp_term = JxW_face[qp] * _tau * (normals[qp] * normals[qp]);
291  // scalar from stabilization term
292  const auto scalar_qp_term = stab_qp_term * scalar_sol[qp];
293  // lm from stabilization term
294  const auto lm_qp_term = stab_qp_term * lm_sol[qp];
295  for (const auto i : index_range(lm_re))
296  lm_re(i) += _lm_phi_face[i][qp] * (scalar_qp_term - vector_qp_term - lm_qp_term);
297  }
298 }
const MooseArray< std::vector< Real > > & _lm_phi_face
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)

◆ makeElemArg()

Moose::ElemArg FunctorInterface::makeElemArg ( const Elem *  elem,
bool  correct_skewnewss = false 
) const
protectedinherited

◆ scalarDirichletJacobian()

void DiffusionLHDGAssemblyHelper::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 
)
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().

384 {
385  for (const auto qp : make_range(qrule_face.n_points()))
386  {
387  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * normals[qp];
388  const auto scalar_qp_term = JxW_face[qp] * _tau * normals[qp] * normals[qp];
389  for (const auto i : index_range(_u_dof_indices))
390  {
391  const auto vector_qpi_term = vector_qp_term * _scalar_phi_face[i][qp];
392  for (const auto j : index_range(_qu_dof_indices))
393  scalar_vector_jac(i, j) -= vector_qpi_term * _vector_phi_face[j][qp];
394 
395  const auto scalar_qpi_term = scalar_qp_term * _scalar_phi_face[i][qp];
396  for (const auto j : index_range(_u_dof_indices))
397  scalar_scalar_jac(i, j) += scalar_qpi_term * _scalar_phi_face[j][qp];
398  }
399  }
400 }
const MooseArray< std::vector< Real > > & _scalar_phi_face
const std::vector< dof_id_type > & _u_dof_indices
const std::vector< dof_id_type > & _qu_dof_indices
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
auto index_range(const T &sizable)

◆ scalarDirichletResidual()

void DiffusionLHDGAssemblyHelper::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 
)
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().

362 {
363  for (const auto qp : make_range(qrule_face.n_points()))
364  {
365  const auto scalar_value = dirichlet_value(
366  Moose::ElemSideQpArg{current_elem, current_side, qp, &qrule_face, q_point_face[qp]},
367  _ti.determineState());
368  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * (vector_sol[qp] * normals[qp]);
369  const auto stab_qp_term = JxW_face[qp] * _tau * normals[qp] * normals[qp];
370  const auto scalar_qp_term = stab_qp_term * scalar_sol[qp];
371  const auto lm_qp_term = stab_qp_term * scalar_value;
372 
373  for (const auto i : index_range(_u_dof_indices))
374  scalar_re(i) += (scalar_qp_term - vector_qp_term - lm_qp_term) * _scalar_phi_face[i][qp];
375  }
376 }
const MooseArray< std::vector< Real > > & _scalar_phi_face
const std::vector< dof_id_type > & _u_dof_indices
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const TransientInterface & _ti
Reference to transient interface.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)
Argument for requesting functor evaluation at quadrature point locations on an element side...

◆ scalarFaceJacobian()

void DiffusionLHDGAssemblyHelper::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 
)
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().

255 {
256  for (const auto qp : make_range(qrule_face.n_points()))
257  {
258  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * normals[qp];
259  const auto stab_qp_term = JxW_face[qp] * _tau * normals[qp] * normals[qp];
260 
261  for (const auto i : make_range(scalar_vector_jac.m()))
262  {
263  const auto vector_qpi_term = vector_qp_term * _scalar_phi_face[i][qp];
264  for (const auto j : make_range(scalar_vector_jac.n()))
265  scalar_vector_jac(i, j) -= vector_qpi_term * _vector_phi_face[j][qp];
266 
267  const auto scalar_qpi_term = stab_qp_term * _scalar_phi_face[i][qp];
268  for (const auto j : make_range(scalar_scalar_jac.n()))
269  scalar_scalar_jac(i, j) += scalar_qpi_term * _scalar_phi_face[j][qp];
270  for (const auto j : make_range(scalar_lm_jac.n()))
271  scalar_lm_jac(i, j) -= scalar_qpi_term * _lm_phi_face[j][qp];
272  }
273  }
274 }
const MooseArray< std::vector< Real > > & _scalar_phi_face
const MooseArray< std::vector< Real > > & _lm_phi_face
unsigned int m() const
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
unsigned int n() const
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face

◆ scalarFaceResidual()

void DiffusionLHDGAssemblyHelper::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 
)
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().

232 {
233  for (const auto qp : make_range(qrule_face.n_points()))
234  {
235  // vector
236  const auto vector_qp_term = JxW_face[qp] * _diff[qp] * (vector_sol[qp] * normals[qp]);
237  // stabilization term
238  const auto stab_qp_term = JxW_face[qp] * _tau * (normals[qp] * normals[qp]);
239  // scalar from stabilization term
240  const auto scalar_qp_term = stab_qp_term * scalar_sol[qp];
241  // lm from stabilization term
242  const auto lm_qp_term = stab_qp_term * lm_sol[qp];
243  for (const auto i : index_range(scalar_re))
244  scalar_re(i) += _scalar_phi_face[i][qp] * (scalar_qp_term - vector_qp_term - lm_qp_term);
245  }
246 }
const MooseArray< std::vector< Real > > & _scalar_phi_face
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
const Real _tau
Our stabilization coefficient.
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)

◆ scalarVolumeJacobian()

void DiffusionLHDGAssemblyHelper::scalarVolumeJacobian ( const MooseArray< Real > &  JxW,
const libMesh::QBase qrule,
DenseMatrix< Number > &  scalar_vector_jac 
)
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().

175 {
176  for (const auto qp : make_range(qrule.n_points()))
177  {
178  const auto qp_term = JxW[qp] * _diff[qp];
179  for (const auto i : make_range(scalar_vector_jac.m()))
180  {
181  const auto qpi_term = qp_term * _grad_scalar_phi[i][qp];
182  // Scalar equation dependence on vector dofs
183  for (const auto j : make_range(scalar_vector_jac.n()))
184  scalar_vector_jac(i, j) += qpi_term * _vector_phi[j][qp];
185  }
186  }
187 }
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
unsigned int m() const
unsigned int n_points() const
const MaterialProperty< Real > & _diff
The diffusivity.
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
unsigned int n() const

◆ scalarVolumeResidual()

void DiffusionLHDGAssemblyHelper::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 
)
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().

152 {
153  for (const auto qp : make_range(qrule.n_points()))
154  {
155  const auto vector_qp_term = JxW[qp] * _diff[qp] * vector_field[qp];
156  // Evaluate source
157  const auto f =
158  source(Moose::ElemQpArg{current_elem, qp, &qrule, q_point[qp]}, _ti.determineState());
159  const auto source_qp_term = JxW[qp] * f;
160 
161  for (const auto i : index_range(scalar_re))
162  {
163  scalar_re(i) += _grad_scalar_phi[i][qp] * vector_qp_term;
164 
165  // Scalar equation RHS
166  scalar_re(i) -= _scalar_phi[i][qp] * source_qp_term;
167  }
168  }
169 }
const MooseArray< std::vector< Real > > & _scalar_phi
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
unsigned int n_points() const
Argument for requesting functor evaluation at a quadrature point location in an element.
const MaterialProperty< Real > & _diff
The diffusivity.
const TransientInterface & _ti
Reference to transient interface.
IntRange< T > make_range(T beg, T end)
auto index_range(const T &sizable)

◆ validParams()

InputParameters DiffusionLHDGAssemblyHelper::validParams ( )
static

Definition at line 22 of file DiffusionLHDGAssemblyHelper.C.

Referenced by DiffusionLHDGDirichletBC::validParams(), DiffusionLHDGKernel::validParams(), and DiffusionLHDGPrescribedGradientBC::validParams().

23 {
24  auto params = emptyInputParameters();
25  params.addRequiredParam<NonlinearVariableName>(
26  "gradient_variable", "The gradient of the diffusing specie concentration");
27  params.addRequiredParam<NonlinearVariableName>(
28  "face_variable", "The concentration of the diffusing specie on faces");
29  params.addRequiredParam<MaterialPropertyName>("diffusivity", "The diffusivity");
30  params.addParam<Real>("tau",
31  1,
32  "The stabilization coefficient required for discontinuous Galerkin "
33  "schemes. This may be set to 0 for a mixed method with Raviart-Thomas.");
34  return params;
35 }
InputParameters emptyInputParameters()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ vectorDirichletResidual()

void DiffusionLHDGAssemblyHelper::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 
)
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().

337 {
338  for (const auto qp : make_range(qrule_face.n_points()))
339  {
340  const auto scalar_value = dirichlet_value(
341  Moose::ElemSideQpArg{current_elem, current_side, qp, &qrule_face, q_point_face[qp]},
342  _ti.determineState());
343  const auto qp_term = JxW_face[qp] * normals[qp] * scalar_value;
344 
345  // External boundary -> Dirichlet faces -> Vector equation RHS
346  for (const auto i : index_range(_qu_dof_indices))
347  vector_re(i) -= qp_term * _vector_phi_face[i][qp];
348  }
349 }
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
const std::vector< dof_id_type > & _qu_dof_indices
unsigned int n_points() const
const TransientInterface & _ti
Reference to transient interface.
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
auto index_range(const T &sizable)
Argument for requesting functor evaluation at quadrature point locations on an element side...

◆ vectorFaceJacobian()

void DiffusionLHDGAssemblyHelper::vectorFaceJacobian ( const MooseArray< Real > &  JxW_face,
const libMesh::QBase qrule_face,
const MooseArray< Point > &  normals,
DenseMatrix< Number > &  vector_lm_jac 
)
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().

210 {
211  for (const auto qp : make_range(qrule_face.n_points()))
212  {
213  const auto qp_term = JxW_face[qp] * normals[qp];
214  // Vector equation dependence on LM dofs
215  for (const auto i : make_range(vector_lm_jac.m()))
216  {
217  const auto qpi_term = qp_term * _vector_phi_face[i][qp];
218  for (const auto j : make_range(vector_lm_jac.n()))
219  vector_lm_jac(i, j) -= qpi_term * _lm_phi_face[j][qp];
220  }
221  }
222 }
const MooseArray< std::vector< Real > > & _lm_phi_face
unsigned int m() const
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
unsigned int n() const
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face

◆ vectorFaceResidual()

void DiffusionLHDGAssemblyHelper::vectorFaceResidual ( const MooseArray< Number > &  lm_sol,
const MooseArray< Real > &  JxW_face,
const libMesh::QBase qrule_face,
const MooseArray< Point > &  normals,
DenseVector< Number > &  vector_re 
)
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().

195 {
196  // Vector equation dependence on LM dofs
197  for (const auto qp : make_range(qrule_face.n_points()))
198  {
199  const auto qp_term = JxW_face[qp] * lm_sol[qp] * normals[qp];
200  for (const auto i : index_range(vector_re))
201  vector_re(i) -= _vector_phi_face[i][qp] * qp_term;
202  }
203 }
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi_face
auto index_range(const T &sizable)

◆ vectorVolumeJacobian()

void DiffusionLHDGAssemblyHelper::vectorVolumeJacobian ( const MooseArray< Real > &  JxW,
const libMesh::QBase qrule,
DenseMatrix< Number > &  vector_vector_jac,
DenseMatrix< Number > &  vector_scalar_jac 
)
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().

128 {
129  for (const auto qp : make_range(qrule.n_points()))
130  for (const auto i : make_range(vector_vector_jac.m()))
131  {
132  // Vector equation dependence on vector dofs
133  const auto vector_qpi_term = JxW[qp] * _vector_phi[i][qp];
134  for (const auto j : make_range(vector_vector_jac.n()))
135  vector_vector_jac(i, j) += vector_qpi_term * _vector_phi[j][qp];
136 
137  // Vector equation dependence on scalar dofs
138  const auto scalar_qpi_term = JxW[qp] * _div_vector_phi[i][qp];
139  for (const auto j : make_range(vector_scalar_jac.n()))
140  vector_scalar_jac(i, j) += scalar_qpi_term * _scalar_phi[j][qp];
141  }
142 }
const MooseArray< std::vector< Real > > & _scalar_phi
const MooseArray< std::vector< Real > > & _div_vector_phi
unsigned int m() const
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
unsigned int n() const

◆ vectorVolumeResidual()

void DiffusionLHDGAssemblyHelper::vectorVolumeResidual ( const MooseArray< Gradient > &  vector_sol,
const MooseArray< Number > &  scalar_sol,
const MooseArray< Real > &  JxW,
const libMesh::QBase qrule,
DenseVector< Number > &  vector_re 
)
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().

107 {
108  for (const auto qp : make_range(qrule.n_points()))
109  {
110  const auto vector_qp_term = JxW[qp] * vector_sol[qp];
111  const auto scalar_qp_term = JxW[qp] * scalar_sol[qp];
112  for (const auto i : index_range(vector_re))
113  {
114  // Vector equation dependence on vector dofs
115  vector_re(i) += _vector_phi[i][qp] * vector_qp_term;
116 
117  // Vector equation dependence on scalar dofs
118  vector_re(i) += _div_vector_phi[i][qp] * scalar_qp_term;
119  }
120  }
121 }
const MooseArray< std::vector< Real > > & _div_vector_phi
unsigned int n_points() const
IntRange< T > make_range(T beg, T end)
const MooseArray< std::vector< RealVectorValue > > & _vector_phi
auto index_range(const T &sizable)

Member Data Documentation

◆ _cached_elem

const Elem* DiffusionLHDGAssemblyHelper::_cached_elem
protected

◆ _dhah_fe_problem

const FEProblemBase& DiffusionLHDGAssemblyHelper::_dhah_fe_problem
private

A reference to the finite element problem used for coupling checks.

Definition at line 281 of file DiffusionLHDGAssemblyHelper.h.

Referenced by checkCoupling().

◆ _dhah_sys

const SystemBase& DiffusionLHDGAssemblyHelper::_dhah_sys
private

A reference to the nonlinear system used for coupling checks.

Definition at line 284 of file DiffusionLHDGAssemblyHelper.h.

Referenced by checkCoupling().

◆ _diff

const MaterialProperty<Real>& DiffusionLHDGAssemblyHelper::_diff
protected

◆ _div_vector_phi

const MooseArray<std::vector<Real> >& DiffusionLHDGAssemblyHelper::_div_vector_phi
protected

Definition at line 249 of file DiffusionLHDGAssemblyHelper.h.

Referenced by vectorVolumeJacobian(), and vectorVolumeResidual().

◆ _grad_scalar_phi

const MooseArray<std::vector<RealVectorValue> >& DiffusionLHDGAssemblyHelper::_grad_scalar_phi
protected

Definition at line 248 of file DiffusionLHDGAssemblyHelper.h.

Referenced by scalarVolumeJacobian(), and scalarVolumeResidual().

◆ _grad_u_var

const MooseVariableFE<RealVectorValue>& DiffusionLHDGAssemblyHelper::_grad_u_var
protected

◆ _lm_lm_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_lm_lm_jac
protected

◆ _lm_phi_face

const MooseArray<std::vector<Real> >& DiffusionLHDGAssemblyHelper::_lm_phi_face
protected

◆ _lm_re

DenseVector<Number> DiffusionLHDGAssemblyHelper::_lm_re
protected

◆ _lm_scalar_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_lm_scalar_jac
protected

◆ _lm_u_dof_indices

const std::vector<dof_id_type>& DiffusionLHDGAssemblyHelper::_lm_u_dof_indices
protected

◆ _lm_u_sol

const MooseArray<Number>& DiffusionLHDGAssemblyHelper::_lm_u_sol
protected

◆ _lm_vector_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_lm_vector_jac
protected

◆ _moose_obj

const MooseObject& DiffusionLHDGAssemblyHelper::_moose_obj
private

A reference to our associated MooseObject for error reporting.

Definition at line 278 of file DiffusionLHDGAssemblyHelper.h.

Referenced by checkCoupling().

◆ _qu_dof_indices

const std::vector<dof_id_type>& DiffusionLHDGAssemblyHelper::_qu_dof_indices
protected

◆ _qu_sol

const MooseArray<libMesh::Gradient>& DiffusionLHDGAssemblyHelper::_qu_sol
protected

◆ _scalar_lm_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_scalar_lm_jac
protected

◆ _scalar_phi

const MooseArray<std::vector<Real> >& DiffusionLHDGAssemblyHelper::_scalar_phi
protected

Definition at line 247 of file DiffusionLHDGAssemblyHelper.h.

Referenced by scalarVolumeResidual(), and vectorVolumeJacobian().

◆ _scalar_phi_face

const MooseArray<std::vector<Real> >& DiffusionLHDGAssemblyHelper::_scalar_phi_face
protected

◆ _scalar_re

DenseVector<Number> DiffusionLHDGAssemblyHelper::_scalar_re
protected

◆ _scalar_scalar_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_scalar_scalar_jac
protected

◆ _scalar_vector_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_scalar_vector_jac
protected

◆ _tau

const Real DiffusionLHDGAssemblyHelper::_tau
protected

◆ _ti

const TransientInterface& DiffusionLHDGAssemblyHelper::_ti
protected

Reference to transient interface.

Definition at line 260 of file DiffusionLHDGAssemblyHelper.h.

Referenced by scalarDirichletResidual(), scalarVolumeResidual(), and vectorDirichletResidual().

◆ _u_dof_indices

const std::vector<dof_id_type>& DiffusionLHDGAssemblyHelper::_u_dof_indices
protected

◆ _u_face_var

const MooseVariableFE<Real>& DiffusionLHDGAssemblyHelper::_u_face_var
protected

◆ _u_sol

const MooseArray<Number>& DiffusionLHDGAssemblyHelper::_u_sol
protected

◆ _u_var

const MooseVariableFE<Real>& DiffusionLHDGAssemblyHelper::_u_var
protected

◆ _vector_lm_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_vector_lm_jac
protected

◆ _vector_phi

const MooseArray<std::vector<RealVectorValue> >& DiffusionLHDGAssemblyHelper::_vector_phi
protected

◆ _vector_phi_face

const MooseArray<std::vector<RealVectorValue> >& DiffusionLHDGAssemblyHelper::_vector_phi_face
protected

◆ _vector_re

DenseVector<Number> DiffusionLHDGAssemblyHelper::_vector_re
protected

◆ _vector_scalar_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_vector_scalar_jac
protected

◆ _vector_vector_jac

DenseMatrix<Number> DiffusionLHDGAssemblyHelper::_vector_vector_jac
protected

The documentation for this class was generated from the following files: