Line data Source code
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 : 15 : /** 16 : * ADInterfaceKernel and ADVectorInterfaceKernel is responsible for interfacing physics across 17 : * subdomains 18 : */ 19 : template <typename T> 20 : class ADInterfaceKernelTempl : public InterfaceKernelBase, public NeighborMooseVariableInterface<T> 21 : { 22 : public: 23 : static InputParameters validParams(); 24 : 25 : ADInterfaceKernelTempl(const InputParameters & parameters); 26 : 27 : /// The primary variable that this interface kernel operates on 28 945 : const MooseVariableFE<T> & variable() const override { return _var; } 29 : 30 : /// The neighbor variable number that this interface kernel operates on 31 821 : const MooseVariableFE<T> & neighborVariable() const override { return _neighbor_var; } 32 : 33 : private: 34 : /** 35 : * Using the passed DGResidual type, selects the correct test function space and residual block, 36 : * and then calls computeQpResidual 37 : */ 38 : void computeElemNeighResidual(Moose::DGResidualType type); 39 : 40 : /** 41 : * Using the passed DGJacobian type, selects the correct test function and trial function spaces 42 : * and 43 : * jacobian block, and then calls computeQpJacobian 44 : */ 45 : void computeElemNeighJacobian(Moose::DGJacobianType type); 46 : 47 : /** 48 : * Using the passed DGJacobian type, selects the correct test function and trial function spaces 49 : * and 50 : * jacobian block, and then calls computeQpOffDiagJacobian with the passed jvar 51 : */ 52 : void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, unsigned int jvar); 53 : 54 : /// Selects the correct Jacobian type and routine to call for the primary variable jacobian 55 : virtual void computeElementOffDiagJacobian(unsigned int jvar) override final; 56 : 57 : /// Selects the correct Jacobian type and routine to call for the secondary variable jacobian 58 : virtual void computeNeighborOffDiagJacobian(unsigned int jvar) override final; 59 : 60 : /// Computes the residual for the current side. 61 : void computeResidual() override final; 62 : 63 : /// Computes the jacobian for the current side. 64 : void computeJacobian() override final; 65 : 66 : protected: 67 : /// Compute residuals at quadrature points 68 : virtual ADReal computeQpResidual(Moose::DGResidualType type) = 0; 69 : 70 : /** 71 : * Put necessary evaluations depending on qp but independent on test functions here 72 : */ 73 2784 : virtual void initQpResidual(Moose::DGResidualType /* type */) {} 74 : 75 : /// The primary side MooseVariable 76 : MooseVariableFE<T> & _var; 77 : 78 : /// Normal vectors at the quadrature points 79 : const MooseArray<Point> & _normals; 80 : 81 : /// Holds the current solution at the current quadrature point on the face. 82 : const ADTemplateVariableValue<T> & _u; 83 : 84 : /// Holds the current solution gradient at the current quadrature point on the face. 85 : const ADTemplateVariableGradient<T> & _grad_u; 86 : 87 : /// The ad version of JxW 88 : const MooseArray<ADReal> & _ad_JxW; 89 : 90 : /// The ad version of coord 91 : const MooseArray<ADReal> & _ad_coord; 92 : 93 : /// The ad version of q_point 94 : const MooseArray<ADPoint> & _ad_q_point; 95 : 96 : /// shape function 97 : const typename OutputTools<T>::VariablePhiValue & _phi; 98 : 99 : /// Side shape function. 100 : const typename OutputTools<T>::VariableTestValue & _test; 101 : 102 : /// Gradient of side shape function 103 : const ADTemplateVariableTestGradient<T> & _grad_test; 104 : 105 : /// Coupled neighbor variable 106 : const MooseVariableFE<T> & _neighbor_var; 107 : 108 : /// Coupled neighbor variable value 109 : const ADTemplateVariableValue<T> & _neighbor_value; 110 : 111 : /// Coupled neighbor variable gradient 112 : const ADTemplateVariableGradient<T> & _grad_neighbor_value; 113 : 114 : /// Side neighbor shape function. 115 : const typename OutputTools<T>::VariablePhiValue & _phi_neighbor; 116 : 117 : /// Side neighbor test function 118 : const typename OutputTools<T>::VariableTestValue & _test_neighbor; 119 : 120 : /// Gradient of side neighbor shape function. No AD because we have not implemented support for 121 : /// neighbor AD FE data in Assembly 122 : const typename OutputTools<T>::VariableTestGradient & _grad_test_neighbor; 123 : 124 : /// Holds residual entries as they are accumulated by this InterfaceKernel 125 : /// This variable is temporarily reserved for RattleSnake 126 : DenseMatrix<Number> _local_kxx; 127 : 128 : private: 129 : Moose::DGResidualType resid_type; 130 : }; 131 : 132 : typedef ADInterfaceKernelTempl<Real> ADInterfaceKernel; 133 : typedef ADInterfaceKernelTempl<RealVectorValue> ADVectorInterfaceKernel;