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 : #include "DGKernelBase.h" 13 : 14 : /** 15 : * The DGKernel class is responsible for calculating the residuals for various 16 : * physics on internal sides (edges/faces). 17 : */ 18 : class DGKernel : public DGKernelBase, public NeighborMooseVariableInterface<Real> 19 : { 20 : public: 21 : /** 22 : * Factory constructor initializes all internal references needed for residual computation. 23 : * 24 : * 25 : * @param parameters The parameters object for holding additional parameters for kernels and 26 : * derived kernels 27 : */ 28 : static InputParameters validParams(); 29 : 30 : DGKernel(const InputParameters & parameters); 31 : 32 : /** 33 : * The variable that this kernel operates on. 34 : */ 35 240109 : virtual const MooseVariableFEBase & variable() const override { return _var; } 36 : 37 : /** 38 : * Computes the residual for this element or the neighbor 39 : */ 40 : virtual void computeElemNeighResidual(Moose::DGResidualType type) override; 41 : 42 : /** 43 : * Computes the element/neighbor-element/neighbor Jacobian 44 : */ 45 : virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override; 46 : 47 : /** 48 : * Computes the element-element off-diagonal Jacobian 49 : */ 50 : virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, 51 : const MooseVariableFEBase & jvar) override; 52 : 53 : protected: 54 : /** 55 : * This is the virtual that derived classes should override for computing the residual on 56 : * neighboring element. 57 : */ 58 : virtual Real computeQpResidual(Moose::DGResidualType type) = 0; 59 : 60 : /** 61 : * This is the virtual that derived classes should override for computing the Jacobian on 62 : * neighboring element. 63 : */ 64 : virtual Real computeQpJacobian(Moose::DGJacobianType type) = 0; 65 : 66 : /** 67 : * This is the virtual that derived classes should override for computing the off-diag Jacobian. 68 : */ 69 : virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar); 70 : 71 : /** 72 : * Insertion point for evaluations that depend on qp but are independent of the test functions. 73 : */ 74 12470944 : virtual void precalculateQpResidual(Moose::DGResidualType /*type*/) {} 75 : 76 : /** 77 : * Insertion point for evaluations that depend on qp but are independent of the test and shape 78 : * functions. 79 : */ 80 1428864 : virtual void precalculateQpJacobian(Moose::DGJacobianType /*type*/) {} 81 : 82 : /** 83 : * Insertion point for evaluations that depend on qp but are independent of the test and shape 84 : * functions for off-diagonal Jacobian assembly. 85 : */ 86 13560 : virtual void precalculateQpOffDiagJacobian(Moose::DGJacobianType /*type*/, 87 : const MooseVariableFEBase & /*jvar*/) 88 : { 89 13560 : } 90 : 91 : /// Variable this kernel operates on 92 : MooseVariable & _var; 93 : /// Holds the current solution at the current quadrature point on the face. 94 : const VariableValue & _u; 95 : /// Holds the current solution gradient at the current quadrature point on the face. 96 : const VariableGradient & _grad_u; 97 : /// Shape functions 98 : const VariablePhiValue & _phi; 99 : /// Gradient of shape function 100 : const VariablePhiGradient & _grad_phi; 101 : /// test functions 102 : const VariableTestValue & _test; 103 : /// Gradient of side shape function 104 : const VariableTestGradient & _grad_test; 105 : /// Side shape function 106 : const VariablePhiValue & _phi_neighbor; 107 : /// Gradient of side shape function 108 : const VariablePhiGradient & _grad_phi_neighbor; 109 : /// Side test function 110 : const VariableTestValue & _test_neighbor; 111 : /// Gradient of side shape function 112 : const VariableTestGradient & _grad_test_neighbor; 113 : /// Holds the current solution at the current quadrature point 114 : const VariableValue & _u_neighbor; 115 : /// Holds the current solution gradient at the current quadrature point 116 : const VariableGradient & _grad_u_neighbor; 117 : };