https://mooseframework.inl.gov
DGKernelBase.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 // local includes
13 #include "MooseArray.h"
14 #include "NeighborResidualObject.h"
15 #include "BlockRestrictable.h"
16 #include "BoundaryRestrictable.h"
19 #include "ElementIDInterface.h"
20 
21 #include "Assembly.h"
22 
23 #define usingDGKernelBaseMembers \
24  usingNeighborCoupleableMembers; \
25  usingFunctionInterfaceMembers; \
26  usingBlockRestrictableMembers; \
27  usingTransientInterfaceMembers
28 
33  public BlockRestrictable,
34  public BoundaryRestrictable,
37  public ElementIDInterface
38 {
39 public:
46 
48 
53 
57  virtual void computeResidual() override;
58 
63 
67  virtual void computeJacobian() override;
68 
73  const MooseVariableFEBase & jvar) = 0;
74 
78  virtual void computeOffDiagJacobian(unsigned int jvar) override;
79 
80  void prepareShapes(unsigned int var_num) override final;
81 
82 protected:
84  const Elem * const & _current_elem;
85 
88 
90  const Elem * const & _neighbor_elem;
91 
94 
96  const unsigned int & _current_side;
98  const Elem *& _current_side_elem;
101 
107  const QBase * const & _qrule;
114 
117  std::vector<MooseVariableFEBase *> _save_in;
118  std::vector<AuxVariableName> _save_in_strings;
119 
122  std::vector<MooseVariableFEBase *> _diag_save_in;
123  std::vector<AuxVariableName> _diag_save_in_strings;
124 
126  const Real & getNeighborElemVolume() const { return _assembly.neighborVolume(); }
127 
128 public:
129  // boundary id used for internal edges (all DG kernels lives on this boundary id -- a made-up
130  // number)
131  static const BoundaryID InternalBndId;
132 
133 protected:
134  unsigned int _qp;
135 
136  unsigned int _i, _j;
137 
139 
142 
144  bool excludeBoundary() const;
145 
146 private:
148  std::set<BoundaryID> _excluded_boundaries;
149 
150  friend class ADDGKernel;
151 };
const Real & neighborVolume()
Returns the reference to the current neighbor volume.
Definition: Assembly.h:469
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
Definition: DGKernelBase.C:160
std::set< BoundaryID > _excluded_boundaries
Broken boundaries.
Definition: DGKernelBase.h:148
const QBase *const & _qrule
Quadrature rule.
Definition: DGKernelBase.h:107
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernelBase.C:26
const Real & _current_side_volume
The volume (or length) of the current side.
Definition: DGKernelBase.h:100
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
const Real & getNeighborElemVolume() const
The volume (or length) of the current neighbor.
Definition: DGKernelBase.h:126
std::vector< AuxVariableName > _save_in_strings
Definition: DGKernelBase.h:118
DGKernelBase(const InputParameters &parameters)
Definition: DGKernelBase.C:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: DGKernelBase.h:121
/class BoundaryRestrictable /brief Provides functionality for limiting the object to certain boundary...
This class provides an interface for common operations on field variables of both FE and FV types wit...
virtual void computeJacobian() override
Computes the jacobian for the current side.
Definition: DGKernelBase.C:139
unsigned int _i
Definition: DGKernelBase.h:136
DGResidualType
Definition: MooseTypes.h:743
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: DGKernelBase.h:122
const MooseArray< Point > & _q_point
Quadrature points.
Definition: DGKernelBase.h:105
virtual void computeElemNeighResidual(Moose::DGResidualType type)=0
Computes the residual for this element or the neighbor.
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
unsigned int _qp
Definition: DGKernelBase.h:134
std::vector< MooseVariableFEBase * > _save_in
Definition: DGKernelBase.h:117
unsigned int _j
Definition: DGKernelBase.h:136
const Moose::CoordinateSystemType & _coord_sys
Coordinate system.
Definition: DGKernelBase.h:103
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
Definition: DGKernelBase.C:209
boundary_id_type BoundaryID
std::vector< AuxVariableName > _diag_save_in_strings
Definition: DGKernelBase.h:123
bool excludeBoundary() const
Check current element if it contains broken boundary.
Definition: DGKernelBase.C:188
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const MooseArray< Real > & _coord
Coordinate transform mainly for curvilinear coordinates.
Definition: DGKernelBase.h:111
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
Definition: DGKernelBase.h:113
BoundaryID _boundary_id
Definition: DGKernelBase.h:138
static const BoundaryID InternalBndId
Definition: DGKernelBase.h:131
static Threads::spin_mutex _resid_vars_mutex
Definition: DGKernelBase.h:140
This is a base class for objects that can provide residual contributions for both local and neighbor ...
const Real & _neighbor_elem_volume
The neighboring element volume.
Definition: DGKernelBase.h:93
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
DGJacobianType
Definition: MooseTypes.h:749
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CoordinateSystemType
Definition: MooseTypes.h:809
virtual void computeResidual() override
Computes the residual for the current side.
Definition: DGKernelBase.C:124
const Elem *const & _current_elem
Current element.
Definition: DGKernelBase.h:84
const MooseArray< Real > & _JxW
Jacobian det times quadrature weighting on quadrature points.
Definition: DGKernelBase.h:109
An interface that restricts an object to subdomains via the &#39;blocks&#39; input parameter.
const unsigned int & _current_side
Current side.
Definition: DGKernelBase.h:96
const Elem *& _current_side_elem
Current side element.
Definition: DGKernelBase.h:98
This interface is designed for DGKernel, InternalSideUserObject, InterfaceUserObject, where material properties on a side of both its primary side (face) and its secondary side (neighbor) all required.
const Elem *const & _neighbor_elem
The neighboring element.
Definition: DGKernelBase.h:90
const InputParameters & parameters() const
Get the parameters of the object.
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)=0
Computes the element/neighbor-element/neighbor Jacobian.
const Real & _current_elem_volume
The volume (or length) of the current element.
Definition: DGKernelBase.h:87
static Threads::spin_mutex _jacoby_vars_mutex
Definition: DGKernelBase.h:141
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: DGKernelBase.h:116
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar)=0
Computes the element-element off-diagonal Jacobian.