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 "IntegratedBCBase.h" 13 : #include "MooseVariableInterface.h" 14 : #include "MooseVariableScalar.h" 15 : 16 : /** 17 : * Base class for deriving any boundary condition of a integrated type 18 : */ 19 : class ArrayIntegratedBC : public IntegratedBCBase, public MooseVariableInterface<RealEigenVector> 20 : { 21 : public: 22 : static InputParameters validParams(); 23 : 24 : ArrayIntegratedBC(const InputParameters & parameters); 25 : 26 8170 : virtual const ArrayMooseVariable & variable() const override { return _var; } 27 : 28 : virtual void computeResidual() override; 29 : virtual void computeJacobian() override; 30 : /** 31 : * Computes d-ivar-residual / d-jvar... 32 : */ 33 : virtual void computeOffDiagJacobian(unsigned int jvar) override; 34 : /** 35 : * Computes jacobian block with respect to a scalar variable 36 : * @param jvar The number of the scalar variable 37 : */ 38 : void computeOffDiagJacobianScalar(unsigned int jvar) override; 39 : 40 : protected: 41 : /** 42 : * Method for computing the residual at quadrature points, to be filled in \p residual. 43 : */ 44 : virtual void computeQpResidual(RealEigenVector & residual) = 0; 45 : 46 : /** 47 : * Method for computing the diagonal Jacobian at quadrature points 48 : */ 49 : virtual RealEigenVector computeQpJacobian(); 50 : 51 : /** 52 : * Method for computing an off-diagonal jacobian component at quadrature points. 53 : */ 54 : virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase & jvar); 55 : 56 : /** 57 : * This is the virtual that derived classes should override for computing a full Jacobian 58 : * component 59 : */ 60 : virtual RealEigenMatrix computeQpOffDiagJacobianScalar(const MooseVariableScalar & jvar); 61 : 62 : /** 63 : * Put necessary evaluations depending on qp but independent on test functions here 64 : */ 65 108400 : virtual void initQpResidual() {} 66 : 67 : /** 68 : * Put necessary evaluations depending on qp but independent on test and shape functions here 69 : */ 70 14464 : virtual void initQpJacobian() {} 71 : 72 : /** 73 : * Put necessary evaluations depending on qp but independent on test and shape functions here for 74 : * off-diagonal Jacobian assembly 75 : */ 76 7552 : virtual void initQpOffDiagJacobian(const MooseVariableFEBase &) {} 77 : 78 : ArrayMooseVariable & _var; 79 : 80 : /// normals at quadrature points 81 : const MooseArray<Point> & _normals; 82 : 83 : /// shape function values (in QPs) 84 : const ArrayVariablePhiValue & _phi; 85 : 86 : /// test function values (in QPs) 87 : const ArrayVariableTestValue & _test; 88 : 89 : /// the values of the unknown variable this BC is acting on 90 : const ArrayVariableValue & _u; 91 : 92 : /// Number of components of the array variable 93 : const unsigned int _count; 94 : 95 : private: 96 : /// Work vector for residual 97 : RealEigenVector _work_vector; 98 : };