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 "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>
25 : InputParameters
26 21736 : ADNodalBCTempl<T, Base>::validParams()
27 : {
28 21736 : return Base::validParams();
29 : }
30 :
31 : template <>
32 : InputParameters
33 3103 : ADNodalBCTempl<RealVectorValue, NodalBCBase>::validParams()
34 : {
35 3103 : InputParameters params = NodalBCBase::validParams();
36 : // The below parameters are useful for vector Nodal BCs
37 12412 : params.addParam<bool>("set_x_comp", true, "Whether to set the x-component of the variable");
38 12412 : params.addParam<bool>("set_y_comp", true, "Whether to set the y-component of the variable");
39 9309 : params.addParam<bool>("set_z_comp", true, "Whether to set the z-component of the variable");
40 3103 : return params;
41 0 : }
42 :
43 : template <>
44 : InputParameters
45 3352 : ADNodalBCTempl<RealVectorValue, ADDirichletBCBase>::validParams()
46 : {
47 3352 : InputParameters params = ADDirichletBCBase::validParams();
48 : // The below parameters are useful for vector Nodal BCs
49 13408 : params.addParam<bool>("set_x_comp", true, "Whether to set the x-component of the variable");
50 13408 : params.addParam<bool>("set_y_comp", true, "Whether to set the y-component of the variable");
51 10056 : params.addParam<bool>("set_z_comp", true, "Whether to set the z-component of the variable");
52 3352 : return params;
53 0 : }
54 :
55 : template <typename T, typename Base>
56 1914 : ADNodalBCTempl<T, Base>::ADNodalBCTempl(const InputParameters & parameters)
57 : : Base(parameters),
58 : MooseVariableInterface<T>(this,
59 : true,
60 : "variable",
61 : Moose::VarKindType::VAR_SOLVER,
62 : std::is_same_v<T, Real> ? Moose::VarFieldType::VAR_FIELD_STANDARD
63 : : std::is_same_v<T, RealVectorValue>
64 : ? Moose::VarFieldType::VAR_FIELD_VECTOR
65 : : Moose::VarFieldType::VAR_FIELD_ARRAY),
66 : ADFunctorInterface(this),
67 3828 : _var(*this->mooseVariable()),
68 1914 : _current_node(_var.node()),
69 3828 : _u(_var.adNodalValue()),
70 7656 : _undisplaced_assembly(_fe_problem.assembly(_tid, _sys.number()))
71 : {
72 : if constexpr (std::is_same_v<T, RealVectorValue>)
73 : {
74 168 : _set_components.resize(Moose::dim);
75 336 : _set_components[0] = this->template getParam<bool>("set_x_comp");
76 336 : _set_components[1] = this->template getParam<bool>("set_y_comp");
77 336 : _set_components[2] = this->template getParam<bool>("set_z_comp");
78 : }
79 : else
80 1746 : _set_components.resize(_var.count(), true);
81 :
82 1914 : _subproblem.haveADObjects(true);
83 :
84 1914 : addMooseVariableDependency(this->mooseVariable());
85 1914 : }
86 :
87 : namespace
88 : {
89 : template <typename T>
90 : const T &
91 333735 : conversionHelper(const T & value, const unsigned int)
92 : {
93 333735 : return value;
94 : }
95 :
96 : template <typename T>
97 : const T &
98 188216 : conversionHelper(const VectorValue<T> & value, const unsigned int i)
99 : {
100 188216 : return value(i);
101 : }
102 :
103 : template <typename T>
104 : const T &
105 27760 : conversionHelper(const Eigen::Matrix<T, Eigen::Dynamic, 1> & value, const unsigned int i)
106 : {
107 27760 : return value(i);
108 : }
109 : }
110 :
111 : template <typename T, typename Base>
112 : void
113 357156 : ADNodalBCTempl<T, Base>::addResidual(const T & residual,
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 803564 : for (const auto i : index_range(dof_indices))
120 446408 : if (_set_components[i])
121 446408 : setResidual(_sys, conversionHelper(residual, i), dof_indices[i]);
122 357156 : }
123 :
124 : template <typename T, typename Base>
125 : template <typename ADResidual>
126 : void
127 84723 : 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 188026 : for (const auto i : index_range(dof_indices))
137 103303 : 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 103303 : addJacobian(_undisplaced_assembly,
141 103303 : Moose::Span(&conversionHelper(residual, i), 1),
142 206606 : Moose::Span(&dof_indices[i], 1),
143 : /*scaling_factor=*/1);
144 84723 : }
145 :
146 : template <typename T, typename Base>
147 : void
148 346396 : ADNodalBCTempl<T, Base>::computeResidual()
149 : {
150 346396 : const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
151 346396 : if (dof_indices.empty())
152 0 : return;
153 :
154 346396 : const auto residual = MetaPhysicL::raw_value(computeQpResidual());
155 :
156 346396 : addResidual(residual, dof_indices);
157 13320 : }
158 :
159 : template <typename T, typename Base>
160 : void
161 73963 : ADNodalBCTempl<T, Base>::computeJacobian()
162 : {
163 73963 : const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
164 73963 : if (dof_indices.empty())
165 0 : return;
166 :
167 73963 : const auto residual = computeQpResidual();
168 :
169 73963 : addJacobian(residual, dof_indices);
170 73963 : }
171 :
172 : template <typename T, typename Base>
173 : void
174 10760 : ADNodalBCTempl<T, Base>::computeResidualAndJacobian()
175 : {
176 10760 : const std::vector<dof_id_type> & dof_indices = _var.dofIndices();
177 10760 : if (dof_indices.empty())
178 0 : return;
179 :
180 10760 : const auto residual = computeQpResidual();
181 :
182 10760 : addResidual(MetaPhysicL::raw_value(residual), dof_indices);
183 10760 : addJacobian(residual, dof_indices);
184 10760 : }
185 :
186 : template <typename T, typename Base>
187 : void
188 76579 : ADNodalBCTempl<T, Base>::computeOffDiagJacobian(const unsigned int jvar_num)
189 : {
190 : // Only need to do this once because AD does all the derivatives at once
191 76579 : if (jvar_num == _var.number())
192 73963 : computeJacobian();
193 76579 : }
194 :
195 : template <typename T, typename Base>
196 : void
197 330 : ADNodalBCTempl<T, Base>::computeOffDiagJacobianScalar(unsigned int)
198 : {
199 : // scalar coupling will have been included in the all-at-once handling in computeOffDiagJacobian
200 330 : }
201 :
202 : template class ADNodalBCTempl<Real, NodalBCBase>;
203 : template class ADNodalBCTempl<RealVectorValue, NodalBCBase>;
204 : template class ADNodalBCTempl<RealEigenVector, NodalBCBase>;
205 : template class ADNodalBCTempl<Real, ADDirichletBCBase>;
206 : template class ADNodalBCTempl<RealVectorValue, ADDirichletBCBase>;
|