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 "KernelBase.h" 13 : #include "MooseVariableInterface.h" 14 : #include "MooseVariableScalar.h" 15 : 16 : class ArrayKernel : public KernelBase, public MooseVariableInterface<RealEigenVector> 17 : { 18 : public: 19 : static InputParameters validParams(); 20 : 21 : ArrayKernel(const InputParameters & parameters); 22 : 23 : /// Compute this ArrayKernel's contribution to the residual 24 : virtual void computeResidual() override; 25 : 26 : /// Compute this ArrayKernel's contribution to the diagonal Jacobian entries 27 : virtual void computeJacobian() override; 28 : 29 : /// Computes full Jacobian of jvar and the array variable this kernel operates on 30 : virtual void computeOffDiagJacobian(unsigned int jvar) override; 31 : 32 : /** 33 : * Computes jacobian block with respect to a scalar variable 34 : * @param jvar The number of the scalar variable 35 : */ 36 : virtual void computeOffDiagJacobianScalar(unsigned int jvar) override; 37 : 38 1110744 : virtual const ArrayMooseVariable & variable() const override { return _var; } 39 : 40 : protected: 41 : /** 42 : * Compute this Kernel's contribution to the residual at the current quadrature point, 43 : * to be filled in \p residual. 44 : */ 45 : virtual void computeQpResidual(RealEigenVector & residual) = 0; 46 : 47 : /** 48 : * Compute this Kernel's contribution to the diagonal Jacobian at the current quadrature point 49 : */ 50 : virtual RealEigenVector computeQpJacobian(); 51 : 52 : /** 53 : * This is the virtual that derived classes should override for computing a full Jacobian 54 : * component 55 : */ 56 : virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase & jvar); 57 : 58 : /** 59 : * This is the virtual that derived classes should override for computing a full Jacobian 60 : * component 61 : */ 62 : virtual RealEigenMatrix computeQpOffDiagJacobianScalar(const MooseVariableScalar & jvar); 63 : 64 : /** 65 : * Put necessary evaluations depending on qp but independent on test functions here 66 : */ 67 13721344 : virtual void initQpResidual() {} 68 : 69 : /** 70 : * Put necessary evaluations depending on qp but independent on test and shape functions here 71 : */ 72 4399538 : virtual void initQpJacobian() {} 73 : 74 : /** 75 : * Put necessary evaluations depending on qp but independent on test and shape functions here for 76 : * off-diagonal Jacobian assembly 77 : */ 78 2640564 : virtual void initQpOffDiagJacobian(const MooseVariableFEBase & jvar) 79 : { 80 2640564 : if (jvar.number() == _var.number()) 81 2559762 : initQpJacobian(); 82 2640564 : } 83 : 84 : /// This is an array kernel so we cast to a ArrayMooseVariable 85 : ArrayMooseVariable & _var; 86 : 87 : /// the current test function 88 : const ArrayVariableTestValue & _test; 89 : 90 : /// gradient of the test function 91 : const ArrayVariableTestGradient & _grad_test; 92 : const MappedArrayVariablePhiGradient & _array_grad_test; 93 : 94 : /// the current shape functions 95 : const ArrayVariablePhiValue & _phi; 96 : 97 : /// gradient of the shape function 98 : const ArrayVariablePhiGradient & _grad_phi; 99 : 100 : /// Holds the solution at current quadrature points 101 : const ArrayVariableValue & _u; 102 : 103 : /// Holds the solution gradient at current quadrature points 104 : const ArrayVariableGradient & _grad_u; 105 : 106 : /// Number of components of the array variable 107 : const unsigned int _count; 108 : 109 : private: 110 : /// Work vector for residual and diag jacobian 111 : RealEigenVector _work_vector; 112 : 113 : /// Work vector for off diag jacobian 114 : RealEigenMatrix _work_matrix; 115 : };