https://mooseframework.inl.gov
ADInterfaceKernel.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 "InterfaceKernelBase.h"
14 
19 template <typename T>
21 {
22 public:
24 
26 
28  const MooseVariableFE<T> & variable() const override { return _var; }
29 
31  const MooseVariableFE<T> & neighborVariable() const override { return _neighbor_var; }
32 
33 private:
39 
46 
53 
55  virtual void computeElementOffDiagJacobian(unsigned int jvar) override final;
56 
58  virtual void computeNeighborOffDiagJacobian(unsigned int jvar) override final;
59 
61  void computeResidual() override final;
62 
64  void computeJacobian() override final;
65 
66 protected:
69 
73  virtual void initQpResidual(Moose::DGResidualType /* type */) {}
74 
77 
80 
83 
86 
89 
92 
95 
98 
101 
104 
107 
110 
113 
116 
119 
123 
127 
128 private:
130 };
131 
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 ADReal computeQpResidual(Moose::DGResidualType type)=0
Compute residuals at quadrature points.
void computeElemNeighJacobian(Moose::DGJacobianType type)
Using the passed DGJacobian type, selects the correct test function and trial function spaces and jac...
ADInterfaceKernelTempl< RealVectorValue > ADVectorInterfaceKernel
const OutputTools< T >::VariablePhiValue & _phi_neighbor
Side neighbor shape function.
const MooseArray< ADPoint > & _ad_q_point
The ad version of q_point.
const MooseArray< ADReal > & _ad_coord
The ad version of coord.
static InputParameters validParams()
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:603
DenseMatrix< Number > _local_kxx
Holds residual entries as they are accumulated by this InterfaceKernel This variable is temporarily r...
const ADTemplateVariableGradient< T > & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void computeResidual() override final
Computes the residual for the current side.
const ADTemplateVariableTestGradient< T > & _grad_test
Gradient of side shape function.
virtual void initQpResidual(Moose::DGResidualType)
Put necessary evaluations depending on qp but independent on test functions here. ...
DGResidualType
Definition: MooseTypes.h:743
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
const ADTemplateVariableValue< T > & _u
Holds the current solution at the current quadrature point on the face.
const MooseVariableFE< T > & _neighbor_var
Coupled neighbor variable.
MooseVariableFE< T > & _var
The primary side MooseVariable.
Enhances MooseVariableInterface interface provide values from neighbor elements.
typename OutputTools< typename Moose::ADType< T >::type >::VariableTestGradient ADTemplateVariableTestGradient
Definition: MooseTypes.h:631
ADInterfaceKernel and ADVectorInterfaceKernel is responsible for interfacing physics across subdomain...
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:606
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
virtual void computeNeighborOffDiagJacobian(unsigned int jvar) override final
Selects the correct Jacobian type and routine to call for the secondary variable jacobian.
const ADTemplateVariableGradient< T > & _grad_neighbor_value
Coupled neighbor variable gradient.
virtual void computeElementOffDiagJacobian(unsigned int jvar) override final
Selects the correct Jacobian type and routine to call for the primary variable jacobian.
DGJacobianType
Definition: MooseTypes.h:749
ADInterfaceKernelTempl< Real > ADInterfaceKernel
void computeElemNeighResidual(Moose::DGResidualType type)
Using the passed DGResidual type, selects the correct test function space and residual block...
const OutputTools< T >::VariableTestValue & _test_neighbor
Side neighbor test function.
const OutputTools< T >::VariableTestGradient & _grad_test_neighbor
Gradient of side neighbor shape function.
const ADTemplateVariableValue< T > & _neighbor_value
Coupled neighbor variable value.
const MooseVariableFE< T > & neighborVariable() const override
The neighbor variable number that this interface kernel operates on.
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
const InputParameters & parameters() const
Get the parameters of the object.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
const OutputTools< T >::VariablePhiValue & _phi
shape function
Moose::DGResidualType resid_type
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const MooseArray< ADReal > & _ad_JxW
The ad version of JxW.
void computeJacobian() override final
Computes the jacobian for the current side.
const OutputTools< T >::VariableTestValue & _test
Side shape function.
const MooseVariableFE< T > & variable() const override
The primary variable that this interface kernel operates on.
ADInterfaceKernelTempl(const InputParameters &parameters)