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 "MooseVariableInterface.h" 13 : #include "NodalBCBase.h" 14 : #include "ADDirichletBCBase.h" 15 : #include "ADFunctorInterface.h" 16 : 17 : /** 18 : * Base class for deriving any automatic differentiation boundary condition of a integrated type 19 : */ 20 : template <typename T, typename Base> 21 : class ADNodalBCTempl : public Base, public MooseVariableInterface<T>, public ADFunctorInterface 22 : { 23 : public: 24 : static InputParameters validParams(); 25 : 26 : ADNodalBCTempl(const InputParameters & parameters); 27 : 28 211675 : const MooseVariableFE<T> & variable() const override { return _var; } 29 : 30 17338 : bool shouldSetComp(unsigned short i) const { return _set_components[i]; } 31 : 32 : protected: 33 : /** 34 : * Compute this NodalBC's contribution to the residual at the current quadrature point 35 : */ 36 : virtual typename Moose::ADType<T>::type computeQpResidual() = 0; 37 : 38 : /// The variable that this NodalBC operates on 39 : MooseVariableFE<T> & _var; 40 : 41 : /// current node being processed 42 : const Node * const & _current_node; 43 : 44 : /// Pseudo-"quadrature point" index (Always zero for the current node) 45 : const unsigned int _qp = 0; 46 : 47 : /// Value of the unknown variable this BC is acting on 48 : const typename Moose::ADType<T>::type & _u; 49 : 50 : const std::array<bool, 3> _set_components; 51 : 52 : using Base::_fe_problem; 53 : using Base::_subproblem; 54 : using Base::_sys; 55 : using Base::_tid; 56 : using Base::addJacobian; 57 : using Base::addMooseVariableDependency; 58 : using Base::setResidual; 59 : 60 : private: 61 : void computeResidual() override final; 62 : void computeJacobian() override final; 63 : void computeResidualAndJacobian() override; 64 : void computeOffDiagJacobian(unsigned int jvar) override final; 65 : void computeOffDiagJacobianScalar(unsigned int jvar) override final; 66 : 67 : /** 68 : * process the residual into the global data structures 69 : */ 70 : template <typename ADResidual> 71 : void addResidual(const ADResidual & residual, const std::vector<dof_id_type> & dof_indices); 72 : 73 : /** 74 : * process the Jacobian into the global data structures 75 : */ 76 : template <typename ADResidual> 77 : void addJacobian(const ADResidual & residual, const std::vector<dof_id_type> & dof_indices); 78 : 79 : /// A reference to the undisplaced assembly in order to ensure data gets correctly incorporated 80 : /// into the global residual/Jacobian 81 : Assembly & _undisplaced_assembly; 82 : }; 83 : 84 : template <> 85 : InputParameters ADNodalBCTempl<RealVectorValue, NodalBCBase>::validParams(); 86 : 87 : template <> 88 : InputParameters ADNodalBCTempl<RealVectorValue, ADDirichletBCBase>::validParams(); 89 : 90 : using ADNodalBC = ADNodalBCTempl<Real, NodalBCBase>; 91 : using ADVectorNodalBC = ADNodalBCTempl<RealVectorValue, NodalBCBase>;