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 : #include "MooseVariableScalar.h" 15 : 16 : /** 17 : * The DGKernel class is responsible for calculating the residuals for various 18 : * physics on internal sides (edges/faces). 19 : */ 20 : class ArrayDGKernel : public DGKernelBase, public NeighborMooseVariableInterface<RealEigenVector> 21 : { 22 : public: 23 : /** 24 : * Factory constructor initializes all internal references needed for residual computation. 25 : * 26 : * 27 : * @param parameters The parameters object for holding additional parameters for kernels and 28 : * derived kernels 29 : */ 30 : static InputParameters validParams(); 31 : 32 : ArrayDGKernel(const InputParameters & parameters); 33 : 34 : /** 35 : * The variable that this kernel operates on. 36 : */ 37 1531 : virtual const MooseVariableFEBase & variable() const override { return _var; } 38 : 39 : /** 40 : * Override this function to consider couplings of components of the array variable 41 : */ 42 : virtual void computeOffDiagJacobian(unsigned int jvar) override; 43 : 44 : /** 45 : * Computes the residual for this element or the neighbor 46 : */ 47 : virtual void computeElemNeighResidual(Moose::DGResidualType type) override; 48 : 49 : /** 50 : * Computes the element/neighbor-element/neighbor Jacobian 51 : */ 52 : virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override; 53 : 54 : /** 55 : * Computes the element-element off-diagonal Jacobian 56 : */ 57 : virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, 58 : const MooseVariableFEBase & jvar) override; 59 : 60 : protected: 61 : /** 62 : * This is the virtual that derived classes should override for computing the residual on 63 : * neighboring element. Residual to be filled in \p residual. 64 : */ 65 : virtual void computeQpResidual(Moose::DGResidualType type, RealEigenVector & residual) = 0; 66 : 67 : /** 68 : * This is the virtual that derived classes should override for computing the Jacobian on 69 : * neighboring element. 70 : */ 71 : virtual RealEigenVector computeQpJacobian(Moose::DGJacobianType); 72 : 73 : /** 74 : * This is the virtual that derived classes should override for computing the off-diag Jacobian. 75 : */ 76 : virtual RealEigenMatrix computeQpOffDiagJacobian(Moose::DGJacobianType type, 77 : const MooseVariableFEBase & jvar); 78 : 79 : /** 80 : * Put necessary evaluations depending on qp but independent on test functions here 81 : */ 82 43296 : virtual void initQpResidual(Moose::DGResidualType) {} 83 : 84 : /** 85 : * Put necessary evaluations depending on qp but independent on test and shape functions here 86 : */ 87 3072 : virtual void initQpJacobian(Moose::DGJacobianType) {} 88 : 89 : /** 90 : * Put necessary evaluations depending on qp but independent on test and shape functions here for 91 : * off-diagonal Jacobian assembly 92 : */ 93 27648 : virtual void initQpOffDiagJacobian(Moose::DGJacobianType, const MooseVariableFEBase &) {} 94 : 95 : /// Variable this kernel operates on 96 : ArrayMooseVariable & _var; 97 : /// Holds the current solution at the current quadrature point on the face. 98 : const ArrayVariableValue & _u; 99 : /// Holds the current solution gradient at the current quadrature point on the face. 100 : const ArrayVariableGradient & _grad_u; 101 : /// Shape functions 102 : const ArrayVariablePhiValue & _phi; 103 : /// Gradient of shape function 104 : const ArrayVariablePhiGradient & _grad_phi; 105 : const MappedArrayVariablePhiGradient & _array_grad_phi; 106 : /// test functions 107 : const ArrayVariableTestValue & _test; 108 : /// Gradient of side shape function 109 : const ArrayVariableTestGradient & _grad_test; 110 : const MappedArrayVariablePhiGradient & _array_grad_test; 111 : /// Side shape function 112 : const ArrayVariablePhiValue & _phi_neighbor; 113 : /// Gradient of side shape function 114 : const ArrayVariablePhiGradient & _grad_phi_neighbor; 115 : const MappedArrayVariablePhiGradient & _array_grad_phi_neighbor; 116 : /// Side test function 117 : const ArrayVariableTestValue & _test_neighbor; 118 : /// Gradient of side shape function 119 : const ArrayVariableTestGradient & _grad_test_neighbor; 120 : const MappedArrayVariablePhiGradient & _array_grad_test_neighbor; 121 : /// Holds the current solution at the current quadrature point 122 : const ArrayVariableValue & _u_neighbor; 123 : /// Holds the current solution gradient at the current quadrature point 124 : const ArrayVariableGradient & _grad_u_neighbor; 125 : /// Normals in type of Eigen map 126 : const std::vector<Eigen::Map<RealDIMValue>> & _array_normals; 127 : /// Number of components of the array variable 128 : const unsigned int _count; 129 : 130 : private: 131 : /// Work vector for residual computation 132 : RealEigenVector _work_vector; 133 : };