https://mooseframework.inl.gov
IPHDGAssemblyHelper.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 
10 #include "IPHDGAssemblyHelper.h"
11 #include "MooseTypes.h"
13 #include "MooseVariableFE.h"
14 #include "MooseVariableScalar.h"
15 #include "SystemBase.h"
16 #include "MooseMesh.h"
17 #include "MooseObject.h"
19 #include "TransientInterface.h"
20 
21 using namespace libMesh;
22 
25 {
26  auto params = emptyInputParameters();
27  params.addRequiredParam<NonlinearVariableName>("face_variable", "The face variable");
28  return params;
29 }
30 
33  const TransientInterface * const ti,
34  SystemBase & sys,
35  const Assembly & assembly,
36  const THREAD_ID tid,
37  const std::set<SubdomainID> & block_ids,
38  const std::set<BoundaryID> & boundary_ids)
39  : ThreeMaterialPropertyInterface(moose_obj, block_ids, boundary_ids),
40  _ti(*ti),
41  _u_var(sys.getFieldVariable<Real>(tid, moose_obj->getParam<NonlinearVariableName>("variable"))),
42  _u_face_var(sys.getFieldVariable<Real>(
43  tid, moose_obj->getParam<NonlinearVariableName>("face_variable"))),
44  _u_dof_indices(_u_var.dofIndices()),
45  _lm_u_dof_indices(_u_face_var.dofIndices()),
46  _u_sol(_u_var.adSln()),
47  _grad_u_sol(_u_var.adGradSln()),
48  _lm_u_sol(_u_face_var.adSln()),
49  _scalar_phi(_u_var.phi()),
50  _grad_scalar_phi(_u_var.gradPhi()),
51  _scalar_phi_face(_u_var.phiFace()),
52  _grad_scalar_phi_face(_u_var.gradPhiFace()),
53  _lm_phi_face(_u_face_var.phiFace()),
54  _elem_volume(assembly.elemVolume()),
55  _side_area(assembly.sideElemVolume()),
56  _ip_current_elem(assembly.elem()),
57  _ip_current_side(assembly.side()),
58  _ip_JxW(assembly.JxW()),
59  _ip_qrule(assembly.qRule()),
60  _ip_q_point(assembly.qPoints()),
61  _ip_JxW_face(assembly.JxWFace()),
62  _ip_qrule_face(assembly.qRuleFace()),
63  _ip_q_point_face(assembly.qPointsFace()),
64  _ip_normals(assembly.normals())
65 {
68 }
69 
70 std::array<ADResidualsPacket, 2>
72 {
75 }
76 
77 std::set<std::string>
79 {
80  return {_u_face_var.name()};
81 }
82 
83 void
85 {
86  for (const auto qp : make_range(_ip_qrule_face->n_points()))
87  {
88  const auto scalar_value = dirichlet_value(
92 
93  for (const auto i : index_range(_lm_re))
94  _lm_re(i) += _ip_JxW_face[qp] * (_lm_u_sol[qp] - scalar_value) * _lm_phi_face[i][qp];
95  }
96 }
97 
98 void
100 {
101  for (const auto qp : make_range(_ip_qrule_face->n_points()))
102  {
103  const auto flux = flux_value(
106  _ti.determineState());
107 
108  for (const auto i : index_range(_lm_re))
109  _lm_re(i) += _ip_JxW_face[qp] * flux * _lm_phi_face[i][qp];
110  }
111 }
DenseVector< ADReal > _scalar_re
Keeps track of stuff related to assembling.
Definition: Assembly.h:101
DenseVector< ADReal > _lm_re
const unsigned int & _ip_current_side
The current element side.
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
void lmPrescribedFlux(const Moose::Functor< Real > &flux_value)
This is a wrapper that forwards calls to the implementation, which can be switched out at any time wi...
const std::string & name() const override
Get the variable name.
const MooseArray< std::vector< Real > > & _lm_phi_face
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void lmDirichlet(const Moose::Functor< Real > &dirichlet_value)
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.
InputParameters emptyInputParameters()
const std::vector< dof_id_type > & _u_dof_indices
const MooseVariableFE< Real > & _u_var
Interface for objects that needs transient capabilities.
const MooseVariableFE< Real > & _u_face_var
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
const MooseArray< Point > & _ip_q_point_face
The physical quadrature point locations on the face.
unsigned int n_points() const
Utility structure for packaging up all of the residual object&#39;s information needed to add into the sy...
const MooseArray< ADReal > & _lm_u_sol
This interface is designed currently for DomainUserObject where material properties on element...
const std::vector< dof_id_type > & _lm_u_dof_indices
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MooseArray< Real > & _ip_JxW_face
The face JxW.
std::set< std::string > additionalROVariables()
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...
std::array< ADResidualsPacket, 2 > taggingData() const
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
unsigned int THREAD_ID
Definition: MooseTypes.h:209
IPHDGAssemblyHelper(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 > &blocks_ids, const std::set< BoundaryID > &boundary_ids)