https://mooseframework.inl.gov
InterfaceKernelBase.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 "BoundaryRestrictable.h"
18 #include "ElementIDInterface.h"
19 
25  public BoundaryRestrictable,
28  public ElementIDInterface
29 {
30 public:
32 
34 
36  virtual const MooseVariableFEBase & neighborVariable() const = 0;
37 
39  virtual void computeElementOffDiagJacobian(unsigned int jvar) = 0;
40 
42  virtual void computeNeighborOffDiagJacobian(unsigned int jvar) = 0;
43 
44  void prepareShapes(unsigned int var_num) override final;
45 
46 protected:
48  virtual Real computeQpJacobian(Moose::DGJacobianType /*type*/) { return 0; }
49 
51  virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType /*type*/, unsigned int /*jvar*/)
52  {
53  return 0;
54  }
55 
57  const Real & getNeighborElemVolume();
58 
60  const Elem * const & _current_elem;
61 
64 
66  const Elem * const & _neighbor_elem;
67 
70 
72  const unsigned int & _current_side;
73 
75  const Elem * const & _current_side_elem;
76 
79 
82 
84  unsigned int _qp;
85 
88 
90  const QBase * const & _qrule;
91 
94 
97 
99  unsigned int _i, _j;
100 
105 
109  std::vector<AuxVariableName> _save_in_strings;
110 
113 
115  std::vector<MooseVariableFEBase *> _primary_save_in_residual_variables;
116 
119 
121  std::vector<MooseVariableFEBase *> _secondary_save_in_residual_variables;
122 
127 
131  std::vector<AuxVariableName> _diag_save_in_strings;
132 
135 
137  std::vector<MooseVariableFEBase *> _primary_save_in_jacobian_variables;
138 
141 
143  std::vector<MooseVariableFEBase *> _secondary_save_in_jacobian_variables;
144 
147 
150 };
static Threads::spin_mutex _resid_vars_mutex
Mutex that prevents multiple threads from saving into the residual aux_var at the same time...
virtual Real computeQpJacobian(Moose::DGJacobianType)
Compute jacobians at quadrature points.
std::vector< MooseVariableFEBase * > _secondary_save_in_residual_variables
The aux variables to save the secondary contributions to.
static InputParameters validParams()
MultiMooseEnum _diag_save_in_var_side
MultiMooseEnum specifying whether jacobian save-in aux variables correspond to primary or secondary s...
const MooseArray< Real > & _JxW
Elemtn Jacobian/quadrature weight.
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
const Real & _neighbor_elem_volume
The neighboring element volume.
const MooseArray< Real > & _coord
Coordinate transformation value; relevant in axisymmetric simulations for example.
const Real & _current_side_volume
The volume (or length) of the current side.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::vector< AuxVariableName > _diag_save_in_strings
The names of the aux variables that will be used to save-in jacobians (includes both primary and seco...
/class BoundaryRestrictable /brief Provides functionality for limiting the object to certain boundary...
virtual const MooseVariableFEBase & neighborVariable() const =0
The neighbor variable number that this interface kernel operates on.
This class provides an interface for common operations on field variables of both FE and FV types wit...
const Elem *const & _current_side_elem
Current side element.
unsigned int _qp
Current quadrature point.
const Real & getNeighborElemVolume()
The volume of the current neighbor.
MultiMooseEnum _save_in_var_side
MultiMooseEnum specifying whether residual save-in aux variables correspond to primary or secondary s...
const Elem *const & _current_elem
Pointer reference to the current element.
bool _has_primary_jacobians_saved_in
Whether there are primary jacobian aux variables.
const MooseArray< Point > & _q_point
Array that holds element quadrature point coordinates.
InterfaceKernelBase(const InputParameters &parameters)
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType, unsigned int)
compute off-diagonal jacobian components at quadrature points
const Elem *const & _neighbor_elem
The neighboring element.
const unsigned int & _current_side
Current side.
std::vector< MooseVariableFEBase * > _primary_save_in_jacobian_variables
The aux variables to save the diagonal Jacobian contributions of the primary variables to...
virtual void computeElementOffDiagJacobian(unsigned int jvar)=0
Selects the correct Jacobian type and routine to call for the primary variable jacobian.
const QBase *const & _qrule
Quadrature rule.
std::vector< MooseVariableFEBase * > _primary_save_in_residual_variables
The aux variables to save the primary residual contributions to.
unsigned int _i
Index for test and trial functions.
const Real & _current_elem_volume
The volume (or length) of the current element.
This is a base class for objects that can provide residual contributions for both local and neighbor ...
virtual void computeNeighborOffDiagJacobian(unsigned int jvar)=0
Selects the correct Jacobian type and routine to call for the secondary variable jacobian.
DGJacobianType
Definition: MooseTypes.h:749
static Threads::spin_mutex _jacoby_vars_mutex
Mutex that prevents multiple threads from saving into the jacobian aux_var at the same time...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CoordinateSystemType
Definition: MooseTypes.h:809
std::vector< MooseVariableFEBase * > _secondary_save_in_jacobian_variables
The aux variables to save the diagonal Jacobian contributions of the secondary variables to...
bool _has_secondary_residuals_saved_in
Whether there are secondary residual aux variables.
std::vector< AuxVariableName > _save_in_strings
The names of the aux variables that will be used to save-in residuals (includes both primary and seco...
bool _has_primary_residuals_saved_in
Whether there are primary residual aux variables.
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 InputParameters & parameters() const
Get the parameters of the object.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
bool _has_secondary_jacobians_saved_in
Whether there are secondary jacobian aux variables.
const Moose::CoordinateSystemType & _coord_sys
Coordinate system.