https://mooseframework.inl.gov
AdvectionIPHDGAssemblyHelper.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "MooseTypes.h"
12 #include "MooseObject.h"
14 
15 using namespace libMesh;
16 
19 {
20  auto params = IPHDGAssemblyHelper::validParams();
21  params.addRequiredParam<MaterialPropertyName>("velocity", "Velocity vector");
22  params.addRequiredParam<Real>("coeff",
23  "A constant coefficient. This could be something like a density");
24  return params;
25 }
26 
28  const MooseObject * const moose_obj,
30  const TransientInterface * const ti,
31  SystemBase & sys,
32  const Assembly & assembly,
33  const THREAD_ID tid,
34  const std::set<SubdomainID> & block_ids,
35  const std::set<BoundaryID> & boundary_ids)
36  : IPHDGAssemblyHelper(moose_obj, mvdi, ti, sys, assembly, tid, block_ids, boundary_ids),
37  _velocity(getADMaterialProperty<RealVectorValue>("velocity")),
38  _face_velocity(getFaceADMaterialProperty<RealVectorValue>("velocity")),
39  _coeff(moose_obj->getParam<Real>("coeff"))
40 {
41 }
42 
43 void
45 {
46  for (const auto qp : make_range(_ip_qrule->n_points()))
47  {
48  const auto adv_quant = _coeff * _u_sol[qp];
49  const auto qp_term = _ip_JxW[qp] * _velocity[qp] * adv_quant;
50  for (const auto i : index_range(_scalar_re))
51  _scalar_re(i) -= _grad_scalar_phi[i][qp] * qp_term;
52  }
53 }
54 
55 ADReal
56 AdvectionIPHDGAssemblyHelper::computeFlux(const unsigned int qp, const ADReal & face_value)
57 {
58  const auto vdotn = _face_velocity[qp] * _ip_normals[qp];
59  const auto face_phi = _coeff * face_value;
60  const auto internal_phi = _coeff * _u_sol[qp];
61  // Short form for writing upwinding. If the velocity is in the direction of the surface normal,
62  // the the internal value is used, else the face value is used
63  return 0.5 * vdotn * (internal_phi + face_phi) + 0.5 * abs(vdotn) * (internal_phi - face_phi);
64 }
65 
66 void
68 {
69  for (const auto qp : make_range(_ip_qrule_face->n_points()))
70  {
71  const auto qp_term = _ip_JxW_face[qp] * computeFlux(qp, _lm_u_sol[qp]);
72  for (const auto i : index_range(_scalar_re))
73  _scalar_re(i) += _scalar_phi_face[i][qp] * qp_term;
74  }
75 }
76 
77 void
79 {
80  for (const auto qp : make_range(_ip_qrule_face->n_points()))
81  {
82  const auto qp_term = _ip_JxW_face[qp] * computeFlux(qp, _lm_u_sol[qp]);
83  for (const auto i : index_range(_lm_re))
84  _lm_re(i) -= _lm_phi_face[i][qp] * qp_term;
85  }
86 }
87 
88 void
90 {
91  for (const auto qp : make_range(_ip_qrule_face->n_points()))
92  {
93  const auto dirichlet_value = dirichlet_functor(
97  const auto qp_term = _ip_JxW_face[qp] * computeFlux(qp, dirichlet_value);
98  for (const auto i : index_range(_scalar_re))
99  _scalar_re(i) += _scalar_phi_face[i][qp] * qp_term;
100  }
101 }
102 
103 void
105 {
106  for (const auto qp : make_range(_ip_qrule_face->n_points()))
107  {
108 #ifndef NDEBUG
109  const auto vdotn = _face_velocity[qp] * _ip_normals[qp];
110  mooseAssert(MetaPhysicL::raw_value(vdotn) >= 0, "The velocity must create outflow conditions");
111 #endif
112  const auto qp_term = _ip_JxW_face[qp] * _coeff * (_lm_u_sol[qp] - _u_sol[qp]);
113  for (const auto i : index_range(_lm_re))
114  // Force the LM solution to be equivalent to the internal solution
115  _lm_re(i) += _lm_phi_face[i][qp] * qp_term;
116  }
117 }
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:50
virtual void scalarDirichlet(const Moose::Functor< Real > &dirichlet_value) override
Weakly imposes a Dirichlet condition for the scalar field in the scalar field equation.
DenseVector< ADReal > _scalar_re
Keeps track of stuff related to assembling.
Definition: Assembly.h:109
ADReal computeFlux(const unsigned int qp, const ADReal &face_value)
compute the face flux, e.g.
DenseVector< ADReal > _lm_re
const unsigned int & _ip_current_side
The current element side.
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
This is a wrapper that forwards calls to the implementation, which can be switched out at any time wi...
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:100
const MooseArray< std::vector< Real > > & _lm_phi_face
void lmOutflow()
prescribes an outflow condition
virtual void lmFace() override
Computes a local residual vector for the weak form: -<Dq*n, > + < * (u - {u}) * n * n...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
Base class for a system (of equations)
Definition: SystemBase.h:85
const QBase *const & _ip_qrule_face
The face qrule.
const TransientInterface & _ti
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:42
const Elem *const & _ip_current_elem
The current element.
const MooseArray< Real > & _ip_JxW
The element JxW.
Base class that declares all the methods for assembling a hybridized interior penalty discontinuous G...
Interface for objects that needs transient capabilities.
virtual void scalarVolume() override
Computes a local residual vector for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity...
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
const MooseArray< Point > & _ip_normals
The normal vector on the face.
const QBase *const & _ip_qrule
The element qrule.
const ADMaterialProperty< RealVectorValue > & _velocity
The velocity in the element interior.
AdvectionIPHDGAssemblyHelper(const MooseObject *const moose_obj, MooseVariableDependencyInterface *const mvdi, const TransientInterface *const ti, SystemBase &sys, const Assembly &assembly, const THREAD_ID tid, const std::set< SubdomainID > &block_ids, const std::set< BoundaryID > &boundary_ids)
const MooseArray< Point > & _ip_q_point_face
The physical quadrature point locations on the face.
unsigned int n_points() const
const ADMaterialProperty< RealVectorValue > & _face_velocity
The velocity on the element faces.
const MooseArray< ADReal > & _lm_u_sol
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void scalarFace() override
Computes a local residual vector for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n...
const MooseArray< ADReal > & _u_sol
const MooseArray< Real > & _ip_JxW_face
The face JxW.
const MooseArray< std::vector< Real > > & _scalar_phi_face
IntRange< T > make_range(T beg, T end)
static InputParameters validParams()
auto index_range(const T &sizable)
Argument for requesting functor evaluation at quadrature point locations on an element side...
const Real _coeff
The advected quantity value is this _coeff value multipled by the variable/side_variable pair (for el...
unsigned int THREAD_ID
Definition: MooseTypes.h:237