www.mooseframework.org
InterfaceKernel.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #ifndef INTERFACEKERNEL_H
16 #define INTERFACEKERNEL_H
17 
18 // local includes
19 #include "MooseArray.h"
20 #include "MooseObject.h"
21 #include "BoundaryRestrictable.h"
22 #include "SetupInterface.h"
23 #include "TransientInterface.h"
24 #include "UserObjectInterface.h"
26 #include "FunctionInterface.h"
27 #include "Restartable.h"
28 #include "ZeroInterface.h"
29 #include "MeshChangedInterface.h"
31 
32 // Forward Declarations
33 class InterfaceKernel;
34 
35 template <>
37 
43  public BoundaryRestrictable,
44  public SetupInterface,
45  public TransientInterface,
46  public FunctionInterface,
47  public UserObjectInterface,
49  public Restartable,
50  public ZeroInterface,
51  public MeshChangedInterface,
53 {
54 public:
56 
58  MooseVariable & variable() const;
59 
61  const MooseVariable & neighborVariable() const;
62 
65 
71 
78 
84  virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, unsigned int jvar);
85 
87  virtual void computeElementOffDiagJacobian(unsigned int jvar);
88 
90  virtual void computeNeighborOffDiagJacobian(unsigned int jvar);
91 
93  virtual void computeResidual();
94 
96  virtual void computeJacobian();
97 
98 protected:
100  virtual Real computeQpResidual(Moose::DGResidualType type) = 0;
101 
103  virtual Real computeQpJacobian(Moose::DGJacobianType type) = 0;
104 
106  virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar);
107 
109  const Real & getNeighborElemVolume();
110 
113 
116 
119 
122 
125 
128 
130  const Elem *& _current_elem;
131 
133  const Real & _current_elem_volume;
134 
136  const Elem *& _neighbor_elem;
137 
139  unsigned int & _current_side;
140 
142  const Elem *& _current_side_elem;
143 
145  const Real & _current_side_volume;
146 
149 
151  unsigned int _qp;
152 
155 
157  QBase *& _qrule;
158 
161 
164 
166  unsigned int _i, _j;
167 
169  const VariableValue & _u;
170 
173 
176 
179 
186 
189 
194 
199 
204 
206  DenseVector<Number> _local_re;
207 
209  DenseMatrix<Number> _local_kxx;
210 
215 
219  std::vector<AuxVariableName> _save_in_strings;
220 
223 
225  std::vector<MooseVariable *> _master_save_in_residual_variables;
226 
229 
231  std::vector<MooseVariable *> _slave_save_in_residual_variables;
232 
237 
241  std::vector<AuxVariableName> _diag_save_in_strings;
242 
245 
247  std::vector<MooseVariable *> _master_save_in_jacobian_variables;
248 
251 
253  std::vector<MooseVariable *> _slave_save_in_jacobian_variables;
254 
256  static Threads::spin_mutex _resid_vars_mutex;
257 
259  static Threads::spin_mutex _jacoby_vars_mutex;
260 };
261 
262 #endif
const VariablePhiGradient & _grad_phi
Shape function gradient.
const VariableTestGradient & _grad_test
Gradient of side shape function.
const VariablePhiValue & _phi
shape function
InterfaceKernel(const InputParameters &parameters)
MultiMooseEnum _save_in_var_side
MultiMooseEnum specifying whether residual save-in aux variables correspond to master or slave side...
const VariablePhiGradient & _grad_phi_neighbor
Gradient of side neighbor shape function.
const Elem *& _neighbor_elem
The neighboring element.
A class for creating restricted objects.
Definition: Restartable.h:31
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, unsigned int jvar)
Using the passed DGJacobian type, selects the correct test function and trial function spaces and jac...
bool _has_slave_residuals_saved_in
Whether there are slave residual aux variables.
Keeps track of stuff related to assembling.
Definition: Assembly.h:63
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
const MooseArray< Real > & _coord
Coordinate transformation value; relevant in axisymmetric simulations for example.
Class for stuff related to variables.
Definition: MooseVariable.h:43
MultiMooseEnum _diag_save_in_var_side
MultiMooseEnum specifying whether jacobian save-in aux variables correspond to master or slave side...
const Elem *& _current_side_elem
Current side element.
const Elem *& _current_elem
Pointer reference to the current element.
virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar)
compute off-diagonal jacobian components at quadrature points
InputParameters validParams< InterfaceKernel >()
const VariableTestGradient & _grad_test_neighbor
Gradient of side neighbor shape function.
MooseVariable & _neighbor_var
Coupled neighbor variable.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this InterfaceKernel.
SubProblem & subProblem()
Return a reference to the subproblem.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _qp
Current quadrature point.
/class BoundaryRestrictable /brief Provides functionality for limiting the object to certain boundary...
virtual Real computeQpResidual(Moose::DGResidualType type)=0
Compute residuals at quadrature points.
const MooseArray< Real > & _JxW
Elemtn Jacobian/quadrature weight.
unsigned int _j
DGResidualType
Definition: MooseTypes.h:184
Base class for a system (of equations)
Definition: SystemBase.h:91
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)
Using the passed DGJacobian type, selects the correct test function and trial function spaces and jac...
std::vector< MooseVariable * > _master_save_in_residual_variables
The aux variables to save the master residual contributions to.
const VariableValue & _neighbor_value
Coupled neighbor variable value.
DenseMatrix< Number > _local_kxx
Holds residual entries as they are accumulated by this InterfaceKernel.
static Threads::spin_mutex _resid_vars_mutex
Mutex that prevents multiple threads from saving into the residual aux_var at the same time...
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
virtual void computeElemNeighResidual(Moose::DGResidualType type)
Using the passed DGResidual type, selects the correct test function space and residual block...
unsigned int _i
Index for test and trial functions.
const VariablePhiValue & _phi_neighbor
Side neighbor shape function.
bool _has_slave_jacobians_saved_in
Whether there are slave jacobian aux variables.
Interface for objects that needs transient capabilities.
const VariableGradient & _grad_neighbor_value
Coupled neighbor variable gradient.
std::vector< AuxVariableName > _diag_save_in_strings
The names of the aux variables that will be used to save-in jacobians (includes both master and slave...
Interface for notifications that the mesh has changed.
const MooseArray< Point > & _q_point
Array that holds element quadrature point coordinates.
virtual void computeResidual()
Computes the residual for the current side.
virtual void computeNeighborOffDiagJacobian(unsigned int jvar)
Selects the correct Jacobian type and routine to call for the slave variable jacobian.
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:36
const Moose::CoordinateSystemType & _coord_sys
Coordinate system.
const VariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
const Real & _current_side_volume
The volume (or length) of the current side.
const Real & getNeighborElemVolume()
The volume of the current neighbor.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:74
SystemBase & _sys
Reference to the nonlinear system.
const MooseVariable & neighborVariable() const
The neighbor variable number that this interface kernel operates on.
Interface for objects that need to use UserObjects.
std::vector< MooseVariable * > _master_save_in_jacobian_variables
The aux variables to save the diagonal Jacobian contributions of the master variables to...
const VariableTestValue & _test
Side shape function.
InterfaceKernel is responsible for interfacing physics across subdomains.
CoordinateSystemType
Definition: MooseTypes.h:212
DGJacobianType
Definition: MooseTypes.h:190
MooseVariable & variable() const
The master variable that this interface kernel operates on.
const VariableTestValue & _test_neighbor
Side neighbor test function.
MooseVariable & _var
The master side MooseVariable.
MatType type
virtual Real computeQpJacobian(Moose::DGJacobianType type)=0
Compute jacobians at quadrature points.
bool _has_master_jacobians_saved_in
Whether there are master jacobian aux variables.
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:53
Interface to bring zero values inside objects.
Definition: ZeroInterface.h:35
Assembly & _assembly
Problem assembly.
MooseMesh & _mesh
The problem mesh.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:53
QBase *& _qrule
Quadrature rule.
THREAD_ID _tid
The thread ID.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
bool _has_master_residuals_saved_in
Whether there are master residual aux variables.
const Real & _current_elem_volume
The volume (or length) of the current element.
const VariableValue & _u
Holds the current solution at the current quadrature point on the face.
static Threads::spin_mutex _jacoby_vars_mutex
Mutex that prevents multiple threads from saving into the jacobian aux_var at the same time...
std::vector< AuxVariableName > _save_in_strings
The names of the aux variables that will be used to save-in residuals (includes both master and slave...
virtual void computeElementOffDiagJacobian(unsigned int jvar)
Selects the correct Jacobian type and routine to call for the master variable jacobian.
std::vector< MooseVariable * > _slave_save_in_jacobian_variables
The aux variables to save the diagonal Jacobian contributions of the slave variables to...
std::vector< MooseVariable * > _slave_save_in_residual_variables
The aux variables to save the slave contributions to.
Interface for objects that need to use functions.
SubProblem & _subproblem
Reference to the controlling finite element problem.
unsigned int & _current_side
Current side.
unsigned int THREAD_ID
Definition: MooseTypes.h:79
virtual void computeJacobian()
Computes the jacobian for the current side.