https://mooseframework.inl.gov
IPHDGAssemblyHelper.h
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 #pragma once
11 
12 #include "Moose.h"
13 #include "MooseTypes.h"
14 #include "MooseArray.h"
15 #include "MooseFunctor.h"
16 #include "Function.h"
17 #include "Kernel.h"
19 #include "TaggingInterface.h"
21 
22 #include "libmesh/vector_value.h"
23 #include <vector>
24 
25 template <typename>
26 class MooseVariableFE;
28 template <typename>
29 class MooseArray;
30 class SystemBase;
31 class MooseMesh;
32 class MooseObject;
34 class TransientInterface;
35 
42 {
43 public:
45 
46  IPHDGAssemblyHelper(const MooseObject * const moose_obj,
48  const TransientInterface * const ti,
49  SystemBase & sys,
50  const Assembly & assembly,
51  const THREAD_ID tid,
52  const std::set<SubdomainID> & blocks_ids,
53  const std::set<BoundaryID> & boundary_ids);
54 
55  void resizeResiduals();
56  virtual void scalarVolume() = 0;
57 
58  //
59  // Methods which can be leveraged both on internal sides in the kernel and by boundary conditions
60  //
61 
62  virtual void scalarFace() = 0;
63 
64  virtual void lmFace() = 0;
65 
66  virtual void scalarDirichlet(const Moose::Functor<Real> & dirichlet_value) = 0;
67 
68  void lmDirichlet(const Moose::Functor<Real> & dirichlet_value);
69 
70  void lmPrescribedFlux(const Moose::Functor<Real> & flux_value);
71 
75  std::array<ADResidualsPacket, 2> taggingData() const;
76 
80  std::set<std::string> additionalROVariables();
81 
82  virtual ~IPHDGAssemblyHelper() = default;
83 
84 protected:
88 
89  // Containers for dof indices
90  const std::vector<dof_id_type> & _u_dof_indices;
91  const std::vector<dof_id_type> & _lm_u_dof_indices;
92 
93  // local solutions at quadrature points
97 
98  // Element shape functions
101 
102  // Face shape functions
106 
109 
111  const Real & _side_area;
112 
114  const Elem * const & _ip_current_elem;
115 
117  const unsigned int & _ip_current_side;
118 
121 
123  const QBase * const & _ip_qrule;
124 
127 
130 
132  const QBase * const & _ip_qrule_face;
133 
136 
139 
140  // Local residual vectors
141  DenseVector<ADReal> _scalar_re, _lm_re;
142 };
143 
144 inline void
146 {
147  _scalar_re.resize(_u_dof_indices.size());
148  _lm_re.resize(_lm_u_dof_indices.size());
149 }
virtual void scalarFace()=0
DenseVector< ADReal > _scalar_re
virtual void lmFace()=0
Keeps track of stuff related to assembling.
Definition: Assembly.h:101
Class for stuff related to variables.
Definition: Adaptivity.h:31
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
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 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)
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
const MooseVariableFE< Real > & _u_var
virtual void scalarVolume()=0
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_normals
The normal vector on the face.
virtual ~IPHDGAssemblyHelper()=default
const QBase *const & _ip_qrule
The element qrule.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
const MooseArray< Point > & _ip_q_point_face
The physical quadrature point locations on the face.
const Real & _elem_volume
The current element volume.
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
forward declarations
const Real & _side_area
The current element side area.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MooseArray< ADReal > & _u_sol
const MooseArray< Real > & _ip_JxW_face
The face JxW.
const MooseArray< std::vector< Real > > & _scalar_phi_face
std::set< std::string > additionalROVariables()
Class for scalar variables (they are different).
const MooseArray< std::vector< RealVectorValue > > & _grad_scalar_phi_face
const MooseArray< Point > & _ip_q_point
The physical quadrature point locations in the element volume.
static InputParameters validParams()
std::array< ADResidualsPacket, 2 > taggingData() const
unsigned int THREAD_ID
Definition: MooseTypes.h:209
virtual void scalarDirichlet(const Moose::Functor< Real > &dirichlet_value)=0
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)
const MooseArray< std::vector< Real > > & _scalar_phi