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 : #include "ADDirichletBCBaseTempl.h" 11 : #include "MooseVariableFE.h" 12 : #include "SystemBase.h" 13 : #include "libmesh/node.h" 14 : 15 : template <typename T> 16 : InputParameters 17 45871 : ADDirichletBCBaseTempl<T>::validParams() 18 : { 19 45871 : InputParameters params = ADNodalBCTempl<T, ADDirichletBCBase>::validParams(); 20 137613 : params.addParam<bool>( 21 91742 : "preset", true, "Whether or not to preset the BC (apply the value before the solve begins)."); 22 45871 : return params; 23 0 : } 24 : 25 : template <typename T> 26 1600 : ADDirichletBCBaseTempl<T>::ADDirichletBCBaseTempl(const InputParameters & parameters) 27 1600 : : ADNodalBCTempl<T, ADDirichletBCBase>(parameters) 28 : { 29 1600 : } 30 : 31 : template <typename T> 32 : void 33 35984 : ADDirichletBCBaseTempl<T>::computeValue(NumericVector<Number> & current_solution) 34 : { 35 : mooseAssert(this->_preset, "BC is not preset"); 36 : 37 35984 : if (_var.isNodalDefined()) 38 : { 39 35984 : const auto n_comp = _current_node->n_comp(_sys.number(), _var.number()); 40 35984 : const auto value = MetaPhysicL::raw_value(computeQpValue()); 41 80637 : for (const auto i : make_range(n_comp)) 42 : { 43 44653 : const auto dof_idx = _current_node->dof_number(_sys.number(), _var.number(), i); 44 : if constexpr (std::is_same<T, Real>::value) 45 : { 46 : mooseAssert(n_comp == 1, "This should only be unity"); 47 27315 : current_solution.set(dof_idx, value); 48 : } 49 : else 50 : { 51 17338 : if (shouldSetComp(i)) 52 17338 : current_solution.set(dof_idx, value(i)); 53 : } 54 : } 55 : } 56 35984 : } 57 : 58 : template <typename T> 59 : typename Moose::ADType<T>::type 60 403552 : ADDirichletBCBaseTempl<T>::computeQpResidual() 61 : { 62 403552 : return _u - computeQpValue(); 63 : } 64 : 65 : template class ADDirichletBCBaseTempl<Real>; 66 : template class ADDirichletBCBaseTempl<RealVectorValue>;