22 #include "libmesh/quadrature.h" 24 template <
typename T,
typename Base>
37 params.
addParam<
bool>(
"set_x_comp",
true,
"Whether to set the x-component of the variable");
38 params.
addParam<
bool>(
"set_y_comp",
true,
"Whether to set the y-component of the variable");
39 params.
addParam<
bool>(
"set_z_comp",
true,
"Whether to set the z-component of the variable");
49 params.
addParam<
bool>(
"set_x_comp",
true,
"Whether to set the x-component of the variable");
50 params.
addParam<
bool>(
"set_y_comp",
true,
"Whether to set the y-component of the variable");
51 params.
addParam<
bool>(
"set_z_comp",
true,
"Whether to set the z-component of the variable");
55 template <
typename T,
typename Base>
67 _var(*this->mooseVariable()),
68 _current_node(_var.node()),
69 _u(_var.adNodalValue()),
70 _undisplaced_assembly(_fe_problem.assembly(_tid, _sys.number()))
72 if constexpr (std::is_same_v<T, RealVectorValue>)
82 _subproblem.haveADObjects(
true);
91 conversionHelper(
const T & value,
const unsigned int)
98 conversionHelper(
const VectorValue<T> & value,
const unsigned int i)
103 template <
typename T>
105 conversionHelper(
const Eigen::Matrix<T, Eigen::Dynamic, 1> & value,
const unsigned int i)
111 template <
typename T,
typename Base>
114 const std::vector<dof_id_type> & dof_indices)
116 mooseAssert(dof_indices.size() <= _set_components.size(),
117 "The number of dof indices must be less than the number of settable components");
120 if (_set_components[i])
121 setResidual(_sys, conversionHelper(residual, i), dof_indices[i]);
124 template <
typename T,
typename Base>
125 template <
typename ADRes
idual>
128 const std::vector<dof_id_type> & dof_indices)
130 mooseAssert(dof_indices.size() <= _set_components.size(),
131 "The number of dof indices must be less than the number of settable components");
132 if (!std::is_same_v<T, RealVectorValue>)
133 mooseAssert(dof_indices.size() == _var.count(),
134 "The number of dof indices should match the variable count");
137 if (_set_components[i])
140 addJacobian(_undisplaced_assembly,
146 template <
typename T,
typename Base>
150 const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
151 if (dof_indices.empty())
156 addResidual(residual, dof_indices);
159 template <
typename T,
typename Base>
163 const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
164 if (dof_indices.empty())
167 const auto residual = computeQpResidual();
169 addJacobian(residual, dof_indices);
172 template <
typename T,
typename Base>
176 const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
177 if (dof_indices.empty())
180 const auto residual = computeQpResidual();
183 addJacobian(residual, dof_indices);
186 template <
typename T,
typename Base>
191 if (jvar_num == _var.number())
195 template <
typename T,
typename Base>
void addJacobian(const ADResidual &residual, const std::vector< dof_id_type > &dof_indices)
process the Jacobian into the global data structures
void computeResidualAndJacobian() override
Base class for deriving any automatic differentiation boundary condition of a integrated type...
void computeJacobian() override final
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
ADNodalBCTempl(const InputParameters ¶meters)
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
InputParameters validParams()
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
MooseVariableFE< T > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
std::vector< bool > _set_components
Replacement for std::span which we only get in c++20.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
VarKindType
Framework-wide stuff.
void computeOffDiagJacobian(unsigned int jvar) override final
static InputParameters validParams()
static InputParameters validParams()
void computeResidual() override final
MooseVariableFE< T > & _var
The variable that this NodalBC operates on.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void computeOffDiagJacobianScalar(unsigned int jvar) override final
Interface for objects that need to get values of MooseVariables.
void addResidual(const T &residual, const std::vector< dof_id_type > &dof_indices)
process the residual into the global data structures
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
auto index_range(const T &sizable)
static InputParameters validParams()