https://mooseframework.inl.gov
ADNodalBC.C
Go to the documentation of this file.
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 "ADNodalBC.h"
11 #include "NodalBCBase.h"
12 #include "ADDirichletBCBase.h"
13 
14 // MOOSE includes
15 #include "Assembly.h"
16 #include "SubProblem.h"
17 #include "SystemBase.h"
18 #include "MooseVariableFE.h"
19 #include "MooseVariableScalar.h"
20 #include "FEProblemBase.h"
21 
22 #include "libmesh/quadrature.h"
23 
24 template <typename T, typename Base>
27 {
28  return Base::validParams();
29 }
30 
31 template <>
34 {
36  // The below parameters are useful for vector Nodal BCs
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");
40  return params;
41 }
42 
43 template <>
46 {
48  // The below parameters are useful for vector Nodal BCs
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");
52  return params;
53 }
54 
55 template <typename T, typename Base>
57  : Base(parameters),
58  MooseVariableInterface<T>(this,
59  true,
60  "variable",
62  std::is_same_v<T, Real> ? Moose::VarFieldType::VAR_FIELD_STANDARD
63  : std::is_same_v<T, RealVectorValue>
66  ADFunctorInterface(this),
67  _var(*this->mooseVariable()),
68  _current_node(_var.node()),
69  _u(_var.adNodalValue()),
70  _undisplaced_assembly(_fe_problem.assembly(_tid, _sys.number()))
71 {
72  if constexpr (std::is_same_v<T, RealVectorValue>)
73  {
75  _set_components[0] = this->template getParam<bool>("set_x_comp");
76  _set_components[1] = this->template getParam<bool>("set_y_comp");
77  _set_components[2] = this->template getParam<bool>("set_z_comp");
78  }
79  else
80  _set_components.resize(_var.count(), true);
81 
82  _subproblem.haveADObjects(true);
83 
84  addMooseVariableDependency(this->mooseVariable());
85 }
86 
87 namespace
88 {
89 template <typename T>
90 const T &
91 conversionHelper(const T & value, const unsigned int)
92 {
93  return value;
94 }
95 
96 template <typename T>
97 const T &
98 conversionHelper(const VectorValue<T> & value, const unsigned int i)
99 {
100  return value(i);
101 }
102 
103 template <typename T>
104 const T &
105 conversionHelper(const Eigen::Matrix<T, Eigen::Dynamic, 1> & value, const unsigned int i)
106 {
107  return value(i);
108 }
109 }
110 
111 template <typename T, typename Base>
112 void
114  const std::vector<dof_id_type> & dof_indices)
115 {
116  mooseAssert(dof_indices.size() <= _set_components.size(),
117  "The number of dof indices must be less than the number of settable components");
118 
119  for (const auto i : index_range(dof_indices))
120  if (_set_components[i])
121  setResidual(_sys, conversionHelper(residual, i), dof_indices[i]);
122 }
123 
124 template <typename T, typename Base>
125 template <typename ADResidual>
126 void
127 ADNodalBCTempl<T, Base>::addJacobian(const ADResidual & residual,
128  const std::vector<dof_id_type> & dof_indices)
129 {
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");
135 
136  for (const auto i : index_range(dof_indices))
137  if (_set_components[i])
138  // If we store into the displaced assembly for nodal bc objects the data never actually makes
139  // it into the global Jacobian
140  addJacobian(_undisplaced_assembly,
141  Moose::Span(&conversionHelper(residual, i), 1),
142  Moose::Span(&dof_indices[i], 1),
143  /*scaling_factor=*/1);
144 }
145 
146 template <typename T, typename Base>
147 void
149 {
150  const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
151  if (dof_indices.empty())
152  return;
153 
154  const auto residual = MetaPhysicL::raw_value(computeQpResidual());
155 
156  addResidual(residual, dof_indices);
157 }
158 
159 template <typename T, typename Base>
160 void
162 {
163  const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
164  if (dof_indices.empty())
165  return;
166 
167  const auto residual = computeQpResidual();
168 
169  addJacobian(residual, dof_indices);
170 }
171 
172 template <typename T, typename Base>
173 void
175 {
176  const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
177  if (dof_indices.empty())
178  return;
179 
180  const auto residual = computeQpResidual();
181 
182  addResidual(MetaPhysicL::raw_value(residual), dof_indices);
183  addJacobian(residual, dof_indices);
184 }
185 
186 template <typename T, typename Base>
187 void
189 {
190  // Only need to do this once because AD does all the derivatives at once
191  if (jvar_num == _var.number())
192  computeJacobian();
193 }
194 
195 template <typename T, typename Base>
196 void
198 {
199  // scalar coupling will have been included in the all-at-once handling in computeOffDiagJacobian
200 }
201 
202 template class ADNodalBCTempl<Real, NodalBCBase>;
VarFieldType
Definition: MooseTypes.h:770
void addJacobian(const ADResidual &residual, const std::vector< dof_id_type > &dof_indices)
process the Jacobian into the global data structures
Definition: ADNodalBC.C:127
void computeResidualAndJacobian() override
Definition: ADNodalBC.C:174
Base class for deriving any automatic differentiation boundary condition of a integrated type...
Definition: ADNodalBC.h:21
void computeJacobian() override final
Definition: ADNodalBC.C:161
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:100
ADNodalBCTempl(const InputParameters &parameters)
Definition: ADNodalBC.C:56
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:162
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
Definition: ADNodalBC.h:50
Replacement for std::span which we only get in c++20.
Definition: MooseTypes.h:1060
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
void computeOffDiagJacobian(unsigned int jvar) override final
Definition: ADNodalBC.C:188
static InputParameters validParams()
Definition: ADNodalBC.C:26
static InputParameters validParams()
void computeResidual() override final
Definition: ADNodalBC.C:148
MooseVariableFE< T > & _var
The variable that this NodalBC operates on.
Definition: ADNodalBC.h:39
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void computeOffDiagJacobianScalar(unsigned int jvar) override final
Definition: ADNodalBC.C:197
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
Definition: ADNodalBC.C:113
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
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()
Definition: NodalBCBase.C:13