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 "ADFunctorInterface.h" 15 : 16 : /** 17 : * Base class for deriving any boundary condition of a integrated type 18 : */ 19 : template <typename T> 20 : class ADIntegratedBCTempl : public IntegratedBCBase, 21 : public MooseVariableInterface<T>, 22 : public ADFunctorInterface 23 : { 24 : public: 25 : static InputParameters validParams(); 26 : 27 : ADIntegratedBCTempl(const InputParameters & parameters); 28 : 29 180746 : const MooseVariableFE<T> & variable() const override { return _var; } 30 : 31 : protected: 32 : void computeResidual() override; 33 : void computeJacobian() override; 34 : void computeResidualAndJacobian() override; 35 : void computeOffDiagJacobian(unsigned int jvar) override; 36 : void computeOffDiagJacobianScalar(unsigned int jvar) override; 37 : 38 : /** 39 : * compute the \p _residuals member for filling the Jacobian. We want to calculate these residuals 40 : * up-front when doing loal derivative indexing because we can use those residuals to fill \p 41 : * _local_ke for every associated jvariable. We do not want to re-do these calculations for every 42 : * jvariable and corresponding \p _local_ke. For global indexing we will simply pass the computed 43 : * \p _residuals directly to \p Assembly::addJacobian 44 : */ 45 : virtual void computeResidualsForJacobian(); 46 : 47 : /** 48 : * Compute this IntegratedBC's contribution to the residual at the current quadrature point 49 : */ 50 : virtual ADReal computeQpResidual() = 0; 51 : 52 : /// The variable that this IntegratedBC operates on 53 : MooseVariableFE<T> & _var; 54 : 55 : /// normals at quadrature points 56 : const MooseArray<ADPoint> & _normals; 57 : 58 : /// (physical) quadrature points 59 : const MooseArray<ADPoint> & _ad_q_points; 60 : 61 : // test functions 62 : 63 : /// test function values (in QPs) 64 : const ADTemplateVariableTestValue<T> & _test; 65 : /// gradients of test functions (in QPs) 66 : const ADTemplateVariableTestGradient<T> & _grad_test; 67 : 68 : /// the values of the unknown variable this BC is acting on 69 : const ADTemplateVariableValue<T> & _u; 70 : /// the gradient of the unknown variable this BC is acting on 71 : const ADTemplateVariableGradient<T> & _grad_u; 72 : 73 : /// The ad version of JxW 74 : const MooseArray<ADReal> & _ad_JxW; 75 : 76 : /// The AD version of coord 77 : const MooseArray<ADReal> & _ad_coord; 78 : 79 : /// The current shape functions 80 : const ADTemplateVariablePhiValue<T> & _phi; 81 : 82 : /// Whether this object is acting on the displaced mesh 83 : const bool _use_displaced_mesh; 84 : 85 : /// Data members for holding residuals 86 : ADReal _r; 87 : std::vector<Real> _residuals; 88 : std::vector<ADReal> _residuals_and_jacobians; 89 : 90 : private: 91 : /** 92 : * compute all the Jacobian entries 93 : */ 94 : void computeADJacobian(); 95 : }; 96 : 97 : using ADIntegratedBC = ADIntegratedBCTempl<Real>; 98 : using ADVectorIntegratedBC = ADIntegratedBCTempl<RealVectorValue>;