https://mooseframework.inl.gov
DiffusionIPHDGAssemblyHelper.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"
13 #include "MooseVariableScalar.h"
14 #include "SystemBase.h"
15 #include "MooseMesh.h"
16 #include "MooseObject.h"
18 
19 using namespace libMesh;
20 
23 {
24  auto params = IPHDGAssemblyHelper::validParams();
25  params.addRequiredParam<MaterialPropertyName>("diffusivity", "The diffusivity");
26  params.addParam<Real>("alpha",
27  1,
28  "The stabilization coefficient required for discontinuous Galerkin "
29  "schemes.");
30  return params;
31 }
32 
34  const MooseObject * const moose_obj,
36  const TransientInterface * const ti,
37  SystemBase & sys,
38  const Assembly & assembly,
39  const THREAD_ID tid,
40  const std::set<SubdomainID> & block_ids,
41  const std::set<BoundaryID> & boundary_ids)
42  : IPHDGAssemblyHelper(moose_obj, mvdi, ti, sys, assembly, tid, block_ids, boundary_ids),
43  _diff(this->getADMaterialProperty<Real>("diffusivity")),
44  _face_diff(this->getFaceADMaterialProperty<Real>("diffusivity")),
45  _alpha(moose_obj->getParam<Real>("alpha"))
46 {
47 }
48 
49 void
51 {
52  for (const auto qp : make_range(_ip_qrule->n_points()))
53  for (const auto i : index_range(_scalar_re))
54  _scalar_re(i) += _ip_JxW[qp] * (_grad_scalar_phi[i][qp] * _diff[qp] * _grad_u_sol[qp]);
55 }
56 
57 void
59 {
60  const auto h_elem = _elem_volume / _side_area;
61 
62  for (const auto i : index_range(_scalar_re))
63  for (const auto qp : make_range(_ip_qrule_face->n_points()))
64  {
65  _scalar_re(i) -= _ip_JxW_face[qp] * _face_diff[qp] * _scalar_phi_face[i][qp] *
66  (_grad_u_sol[qp] * _ip_normals[qp]);
67  _scalar_re(i) -= _ip_JxW_face[qp] * _alpha / h_elem * _face_diff[qp] *
68  (_lm_u_sol[qp] - _u_sol[qp]) * _scalar_phi_face[i][qp];
69  _scalar_re(i) += _ip_JxW_face[qp] * (_lm_u_sol[qp] - _u_sol[qp]) * _face_diff[qp] *
71  }
72 }
73 
74 void
76 {
77  const auto h_elem = _elem_volume / _side_area;
78 
79  for (const auto i : index_range(_lm_re))
80  for (const auto qp : make_range(_ip_qrule_face->n_points()))
81  {
82  _lm_re(i) += _ip_JxW_face[qp] * _face_diff[qp] * _grad_u_sol[qp] * _ip_normals[qp] *
83  _lm_phi_face[i][qp];
84  _lm_re(i) += _ip_JxW_face[qp] * _alpha / h_elem * _face_diff[qp] *
85  (_lm_u_sol[qp] - _u_sol[qp]) * _lm_phi_face[i][qp];
86  }
87 }
88 
89 void
91 {
92  const auto h_elem = _elem_volume / _side_area;
93 
94  for (const auto qp : make_range(_ip_qrule_face->n_points()))
95  {
96  const auto scalar_value = dirichlet_value(
100 
101  for (const auto i : index_range(_u_dof_indices))
102  {
103  _scalar_re(i) -= _ip_JxW_face[qp] * _face_diff[qp] * _scalar_phi_face[i][qp] *
104  (_grad_u_sol[qp] * _ip_normals[qp]);
105  _scalar_re(i) -= _ip_JxW_face[qp] * _alpha / h_elem * _face_diff[qp] *
106  (scalar_value - _u_sol[qp]) * _scalar_phi_face[i][qp];
107  _scalar_re(i) += _ip_JxW_face[qp] * (scalar_value - _u_sol[qp]) * _face_diff[qp] *
108  _grad_scalar_phi_face[i][qp] * _ip_normals[qp];
109  }
110  }
111 }
const ADMaterialProperty< Real > & _diff
The diffusivity in the element volume.
DenseVector< ADReal > _scalar_re
Keeps track of stuff related to assembling.
Definition: Assembly.h:101
const Real _alpha
Our stabilization coefficient.
virtual void scalarVolume() override
Computes a local residual vector for the weak form: (Dq, grad(w)) - (f, w) where D is the diffusivity...
DenseVector< ADReal > _lm_re
const MooseArray< ADRealVectorValue > & _grad_u_sol
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...
const MooseArray< std::vector< Real > > & _lm_phi_face
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:84
const QBase *const & _ip_qrule_face
The face qrule.
const TransientInterface & _ti
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...
const std::vector< dof_id_type > & _u_dof_indices
Interface for objects that needs transient capabilities.
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 MooseArray< Point > & _ip_q_point_face
The physical quadrature point locations on the face.
const Real & _elem_volume
The current element volume.
unsigned int n_points() const
const ADMaterialProperty< Real > & _face_diff
The diffusivity on the element faces.
const MooseArray< ADReal > & _lm_u_sol
const Real & _side_area
The current element side area.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void lmFace() override
Computes a local residual vector for the weak form: -<Dq*n, > + < * (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)
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi_face
DiffusionIPHDGAssemblyHelper(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)
virtual void scalarFace() override
Computes a local residual vector for the weak form: -<Dq*n, w> + < * (u - {u}) * n * n...
virtual void scalarDirichlet(const Moose::Functor< Real > &dirichlet_value) override
Weakly imposes a Dirichlet condition for the scalar field in the scalar field equation.
static InputParameters validParams()
auto index_range(const T &sizable)
Argument for requesting functor evaluation at quadrature point locations on an element side...
unsigned int THREAD_ID
Definition: MooseTypes.h:209