www.mooseframework.org
InterfaceKernel.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "InterfaceKernelBase.h"
14 
15 #define TemplateVariableValue typename OutputTools<T>::VariableValue
16 #define TemplateVariableGradient typename OutputTools<T>::VariableGradient
17 #define TemplateVariablePhiValue typename OutputTools<T>::VariablePhiValue
18 #define TemplateVariablePhiGradient typename OutputTools<T>::VariablePhiGradient
19 #define TemplateVariableTestValue typename OutputTools<T>::VariableTestValue
20 #define TemplateVariableTestGradient typename OutputTools<T>::VariableTestGradient
21 
22 // Forward Declarations
23 template <typename T>
25 
28 
33 template <typename T>
35 {
36 public:
38 
40 
42  virtual const MooseVariableFE<T> & variable() const override { return _var; }
43 
45  virtual const MooseVariableFE<T> & neighborVariable() const override { return _neighbor_var; }
46 
52 
59 
65  virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, unsigned int jvar);
66 
68  virtual void computeElementOffDiagJacobian(unsigned int jvar) override;
69 
71  virtual void computeNeighborOffDiagJacobian(unsigned int jvar) override;
72 
74  virtual void computeResidual() override;
75 
77  virtual void computeJacobian() override;
78 
80  virtual void computeResidualAndJacobian() override;
81 
84 
88  virtual void initQpResidual(Moose::DGResidualType /* type */) {}
89 
93  virtual void initQpJacobian(Moose::DGJacobianType /* type */) {}
94 
99  virtual void initQpOffDiagJacobian(Moose::DGJacobianType /* type */, unsigned int /* jvar */) {}
100 
101 protected:
104 
107 
109  const TemplateVariableValue & _u;
110 
112  const TemplateVariableGradient & _grad_u;
113 
115  const TemplateVariablePhiValue & _phi;
116 
118  const TemplateVariablePhiGradient & _grad_phi;
119 
121  const TemplateVariableTestValue & _test;
122 
124  const TemplateVariableTestGradient & _grad_test;
125 
128 
130  const TemplateVariableValue & _neighbor_value;
131 
133  const TemplateVariableGradient & _grad_neighbor_value;
134 
136  const TemplateVariablePhiValue & _phi_neighbor;
137 
139  const TemplateVariablePhiGradient & _grad_phi_neighbor;
140 
142  const TemplateVariableTestValue & _test_neighbor;
143 
145  const TemplateVariableTestGradient & _grad_test_neighbor;
146 
150 };
virtual void initQpOffDiagJacobian(Moose::DGJacobianType, unsigned int)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
const TemplateVariableValue & _neighbor_value
Coupled neighbor variable value.
const TemplateVariableValue & _u
Holds the current solution at the current quadrature point on the face.
InterfaceKernel and VectorInterfaceKernel is responsible for interfacing physics across subdomains...
InterfaceKernelTempl(const InputParameters &parameters)
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const TemplateVariableTestGradient & _grad_test
Gradient of side shape function.
const TemplateVariablePhiValue & _phi
shape function
DGResidualType
Definition: MooseTypes.h:656
const TemplateVariableTestGradient & _grad_test_neighbor
Gradient of side neighbor shape function.
virtual void initQpResidual(Moose::DGResidualType)
Put necessary evaluations depending on qp but independent on test functions here. ...
Enhances MooseVariableInterface interface provide values from neighbor elements.
virtual const MooseVariableFE< T > & variable() const override
The primary variable that this interface kernel operates on.
const MooseVariableFE< T > & _neighbor_var
Coupled neighbor variable.
virtual void computeElementOffDiagJacobian(unsigned int jvar) override
Selects the correct Jacobian type and routine to call for the primary variable jacobian.
virtual void computeElemNeighResidual(Moose::DGResidualType type)
Using the passed DGResidual type, selects the correct test function space and residual block...
virtual void computeResidualAndJacobian() override
Computes the residual and Jacobian for the current side.
virtual const MooseVariableFE< T > & neighborVariable() const override
The neighbor variable number that this interface kernel operates on.
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...
virtual void computeJacobian() override
Computes the jacobian for the current side.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const TemplateVariableTestValue & _test
Side shape function.
const TemplateVariablePhiGradient & _grad_phi_neighbor
Gradient of side neighbor shape function.
DGJacobianType
Definition: MooseTypes.h:662
virtual void computeElemNeighJacobian(Moose::DGJacobianType type)
Using the passed DGJacobian type, selects the correct test function and trial function spaces and jac...
InterfaceKernelTempl< RealVectorValue > VectorInterfaceKernel
const TemplateVariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void initQpJacobian(Moose::DGJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
DenseMatrix< Number > _local_kxx
Holds residual entries as they are accumulated by this InterfaceKernel This variable is temporarily r...
InterfaceKernelTempl< Real > InterfaceKernel
const InputParameters & parameters() const
Get the parameters of the object.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
const TemplateVariablePhiGradient & _grad_phi
Shape function gradient.
MooseVariableFE< T > & _var
The primary side MooseVariable.
virtual void computeResidual() override
Computes the residual for the current side.
const TemplateVariableTestValue & _test_neighbor
Side neighbor test function.
const TemplateVariablePhiValue & _phi_neighbor
Side neighbor shape function.
static InputParameters validParams()
virtual void computeNeighborOffDiagJacobian(unsigned int jvar) override
Selects the correct Jacobian type and routine to call for the secondary variable jacobian.
const TemplateVariableGradient & _grad_neighbor_value
Coupled neighbor variable gradient.
virtual Real computeQpResidual(Moose::DGResidualType type)=0
Compute residuals at quadrature points.