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 "Coupleable.h"
11 : #include "Problem.h"
12 : #include "SubProblem.h"
13 : #include "FEProblem.h"
14 : #include "MooseVariableScalar.h"
15 : #include "MooseVariableFE.h"
16 : #include "InputParameters.h"
17 : #include "MooseObject.h"
18 : #include "SystemBase.h"
19 : #include "AuxiliarySystem.h"
20 :
21 : #include "AuxKernel.h"
22 : #include "ElementUserObject.h"
23 : #include "NodalUserObject.h"
24 : #include "NodeFaceConstraint.h"
25 : #include "NodeElemConstraintBase.h"
26 :
27 385799 : Coupleable::Coupleable(const MooseObject * moose_object, bool nodal, bool is_fv)
28 771598 : : _c_parameters(moose_object->parameters()),
29 385799 : _c_name(moose_object->name()),
30 385799 : _c_type(moose_object->type()),
31 1543196 : _c_fe_problem(*_c_parameters.getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
32 771598 : _c_sys(_c_parameters.isParamValid("_sys") ? _c_parameters.get<SystemBase *>("_sys") : nullptr),
33 385799 : _new_to_deprecated_coupled_vars(_c_parameters.getNewToDeprecatedVarMap()),
34 385799 : _c_nodal(nodal),
35 385799 : _c_is_implicit(_c_parameters.have_parameter<bool>("implicit")
36 385799 : ? _c_parameters.get<bool>("implicit")
37 : : true),
38 385799 : _c_allow_element_to_nodal_coupling(
39 385799 : _c_parameters.have_parameter<bool>("_allow_nodal_to_elemental_coupling")
40 385799 : ? _c_parameters.get<bool>("_allow_nodal_to_elemental_coupling")
41 : : false),
42 385799 : _c_tid(_c_parameters.get<THREAD_ID>("_tid")),
43 385799 : _zero(_c_fe_problem._zero[_c_tid]),
44 385799 : _phi_zero(_c_fe_problem._phi_zero[_c_tid]),
45 385799 : _ad_zero(_c_fe_problem._ad_zero[_c_tid]),
46 385799 : _grad_zero(_c_fe_problem._grad_zero[_c_tid]),
47 385799 : _ad_grad_zero(_c_fe_problem._ad_grad_zero[_c_tid]),
48 385799 : _grad_phi_zero(_c_fe_problem._grad_phi_zero[_c_tid]),
49 385799 : _second_zero(_c_fe_problem._second_zero[_c_tid]),
50 385799 : _ad_second_zero(_c_fe_problem._ad_second_zero[_c_tid]),
51 385799 : _second_phi_zero(_c_fe_problem._second_phi_zero[_c_tid]),
52 385799 : _vector_zero(_c_fe_problem._vector_zero[_c_tid]),
53 385799 : _vector_curl_zero(_c_fe_problem._vector_curl_zero[_c_tid]),
54 385799 : _coupleable_neighbor(_c_parameters.have_parameter<bool>("_neighbor")
55 385799 : ? _c_parameters.get<bool>("_neighbor")
56 : : false),
57 385799 : _coupleable_max_qps(Moose::constMaxQpsPerElem),
58 385799 : _is_fv(is_fv),
59 385799 : _obj(moose_object),
60 3857990 : _writable_coupled_variables(libMesh::n_threads())
61 : {
62 1543196 : SubProblem & problem = *_c_parameters.getCheckedPointerParam<SubProblem *>("_subproblem");
63 385799 : _obj->getMooseApp().registerInterfaceObject(*this);
64 :
65 385799 : unsigned int optional_var_index_counter = 0;
66 :
67 : // Coupling
68 674873 : for (auto iter = _c_parameters.coupledVarsBegin(); iter != _c_parameters.coupledVarsEnd(); ++iter)
69 : {
70 289082 : std::string name = *iter;
71 :
72 289082 : std::vector<std::string> vars = _c_parameters.getVecMooseType(name);
73 289082 : if (vars.size() > 0)
74 : {
75 161983 : for (const auto & coupled_var_name : vars)
76 : {
77 88713 : if (problem.hasVariable(coupled_var_name))
78 : {
79 : MooseVariableFieldBase * moose_var =
80 87122 : &problem.getVariable(_c_tid,
81 : coupled_var_name,
82 : Moose::VarKindType::VAR_ANY,
83 87122 : Moose::VarFieldType::VAR_FIELD_ANY);
84 87122 : _coupled_vars[name].push_back(moose_var);
85 87122 : _coupled_moose_vars.push_back(moose_var);
86 87122 : if (auto * tmp_var = dynamic_cast<MooseVariable *>(moose_var))
87 80015 : _coupled_standard_moose_vars.push_back(tmp_var);
88 7107 : else if (auto * tmp_var = dynamic_cast<VectorMooseVariable *>(moose_var))
89 1118 : _coupled_vector_moose_vars.push_back(tmp_var);
90 5989 : else if (auto * tmp_var = dynamic_cast<ArrayMooseVariable *>(moose_var))
91 1579 : _coupled_array_moose_vars.push_back(tmp_var);
92 4410 : else if (auto * tmp_var = dynamic_cast<MooseVariableFV<Real> *>(moose_var))
93 : {
94 : // We are using a finite volume variable through add*CoupledVar as opposed to getFunctor
95 : // so we can be reasonably confident that the variable values will be obtained using
96 : // traditional pre-evaluation and quadrature point indexing
97 4056 : tmp_var->requireQpComputations();
98 4056 : _coupled_fv_moose_vars.push_back(tmp_var);
99 : }
100 354 : else if (auto * tmp_var = dynamic_cast<MooseLinearVariableFV<Real> *>(moose_var))
101 354 : _coupled_fv_moose_vars.push_back(tmp_var);
102 : else
103 0 : _obj->paramError(name, "provided c++ type for variable parameter is not supported");
104 : }
105 1591 : else if (problem.hasScalarVariable(coupled_var_name))
106 : {
107 : MooseVariableScalar * moose_scalar_var =
108 1583 : &problem.getScalarVariable(_c_tid, coupled_var_name);
109 1583 : _c_coupled_scalar_vars[name].push_back(moose_scalar_var);
110 : }
111 : else
112 8 : _obj->paramError(name, "coupled variable '", coupled_var_name, "' was not found");
113 : }
114 : }
115 : else // This means it was optional coupling. Let's assign a unique id to this variable
116 : {
117 215804 : _optional_var_index[name].assign(_c_parameters.numberDefaultCoupledValues(name), 0);
118 432311 : for (unsigned int j = 0; j < _optional_var_index[name].size(); ++j)
119 216507 : _optional_var_index[name][j] =
120 216507 : std::numeric_limits<unsigned int>::max() - optional_var_index_counter;
121 215804 : ++optional_var_index_counter;
122 : }
123 289074 : }
124 771590 : }
125 :
126 : bool
127 297265 : Coupleable::isCoupled(const std::string & var_name_in, unsigned int i) const
128 : {
129 297265 : const auto var_name = _c_parameters.checkForRename(var_name_in);
130 :
131 297265 : auto it = _coupled_vars.find(var_name);
132 297265 : if (it != _coupled_vars.end())
133 166015 : return (i < it->second.size());
134 : else
135 : {
136 : // Make sure the user originally requested this value in the InputParameter syntax
137 131250 : if (!_c_parameters.hasCoupledValue(var_name))
138 4 : mooseError(_c_name,
139 : ": The coupled variable \"",
140 : var_name,
141 : "\" was never added to this object's "
142 : "InputParameters, please double-check your "
143 : "spelling");
144 :
145 131246 : return false;
146 : }
147 297261 : }
148 :
149 : bool
150 1616 : Coupleable::isCoupledConstant(const std::string & var_name) const
151 : {
152 1616 : return _c_parameters.hasDefaultCoupledValue(var_name);
153 : }
154 :
155 : unsigned int
156 97121 : Coupleable::coupledComponents(const std::string & var_name_in) const
157 : {
158 97121 : const auto var_name = _c_parameters.checkForRename(var_name_in);
159 :
160 97121 : if (isCoupled(var_name))
161 : {
162 : mooseAssert(_coupled_vars.find(var_name) != _coupled_vars.end(),
163 : var_name << " must not actually be coupled!");
164 4504 : return _coupled_vars.at(var_name).size();
165 : }
166 : else
167 : {
168 92617 : if (_c_parameters.hasDefaultCoupledValue(var_name))
169 72 : return _c_parameters.numberDefaultCoupledValues(var_name);
170 : else
171 92545 : return 0;
172 : }
173 97121 : }
174 :
175 : void
176 162811 : checkComponent(const MooseObject * obj,
177 : unsigned int comp,
178 : unsigned int bound,
179 : const std::string & var_name)
180 : {
181 162811 : if (bound > 0 && comp >= bound)
182 8 : obj->paramError(
183 : var_name, "component ", comp, " is out of range for this variable (max ", bound - 1, ")");
184 162803 : }
185 :
186 : // calls to this must go *after* get[bla]Var calls and (checking for nullptr
187 : // return). Because checkFuncType calls coupledCallback which should only be
188 : // called if the variables was actually coupled.
189 : void
190 134337 : Coupleable::checkFuncType(const std::string var_name, VarType t, FuncAge age) const
191 : {
192 134337 : if (t == VarType::Gradient && _c_nodal)
193 0 : mooseError(_c_name, ": nodal variables do not have gradients at nodes");
194 :
195 134337 : if (age == FuncAge::Old || age == FuncAge::Older || t == VarType::GradientDot ||
196 : t == VarType::Dot)
197 51366 : validateExecutionerType(var_name, "coupled[Vector][Gradient/Dot]Old[er]");
198 134329 : if (age == FuncAge::Older && !_c_is_implicit)
199 0 : mooseError("object '",
200 0 : _c_name,
201 : "' uses older variable values that are unavailable with explicit schemes");
202 :
203 134329 : coupledCallback(var_name, age == FuncAge::Old || age == FuncAge::Older);
204 134329 : }
205 :
206 : bool
207 163417 : Coupleable::checkVar(const std::string & var_name_in,
208 : unsigned int comp,
209 : unsigned int comp_bound) const
210 : {
211 163417 : const auto var_name = _c_parameters.checkForRename(var_name_in);
212 163417 : auto it = _c_coupled_scalar_vars.find(var_name);
213 163417 : if (it != _c_coupled_scalar_vars.end())
214 : {
215 4 : std::string cvars;
216 8 : for (auto jt : it->second)
217 4 : cvars += " " + jt->name();
218 :
219 4 : _obj->paramError(var_name,
220 : "cannot couple '",
221 : var_name,
222 : "' to a scalar variable (",
223 : cvars,
224 : ") where field variable is expected");
225 0 : }
226 :
227 163413 : if (!isCoupled(var_name, comp))
228 2060 : return false; // return false since variable is *not* coupled
229 :
230 161349 : auto vars_vector_it = _coupled_vars.find(var_name);
231 161349 : if (vars_vector_it == _coupled_vars.end())
232 0 : mooseError(_c_name, ": Trying to get a coupled var ", var_name, " that doesn't exist");
233 :
234 161349 : const auto & vars_vector = vars_vector_it->second;
235 :
236 161349 : auto bound = comp_bound ? comp_bound : vars_vector.size();
237 161349 : checkComponent(_obj, comp, bound, var_name);
238 :
239 : // We should know we have a variable now
240 161349 : const auto * var = vars_vector[comp];
241 161349 : if (!var)
242 0 : mooseError(
243 0 : _c_name,
244 : ": We did all our checks for the existence of a var, yet we still don't have a var!?");
245 :
246 : // Only perform the following checks for objects that feed into residuals/Jacobians, e.g. objects
247 : // that inherit from the TaggingInterface
248 161349 : if (_c_parameters.have_parameter<MultiMooseEnum>("vector_tags"))
249 : {
250 : // Are we attempting to couple to a non-FV var in an FV object?
251 18045 : if (!var->isFV() && _is_fv)
252 0 : mooseError("Attempting to couple non-FV variable ",
253 0 : var->name(),
254 : " into an FV object ",
255 0 : _c_name,
256 : ". This is not currently supported");
257 : }
258 :
259 161349 : if (!(vars_vector[comp])->isNodal() && _c_nodal && !_c_allow_element_to_nodal_coupling)
260 8 : mooseError(_c_name, ": cannot couple elemental variables into nodal objects");
261 :
262 161341 : return true;
263 163401 : }
264 :
265 : const MooseVariableFieldBase *
266 0 : Coupleable::getFEVar(const std::string & var_name, unsigned int comp) const
267 : {
268 0 : mooseDeprecated("Coupleable::getFEVar is deprecated. Please use Coupleable::getFieldVar instead. "
269 : "Note that this method could potentially return a finite volume variable");
270 0 : return getFieldVar(var_name, comp);
271 : }
272 :
273 : MooseVariableFieldBase *
274 7154 : Coupleable::getFieldVar(const std::string & var_name, unsigned int comp)
275 : {
276 7154 : return getVarHelper<MooseVariableFieldBase>(var_name, comp);
277 : }
278 :
279 : const MooseVariableFieldBase *
280 23314 : Coupleable::getFieldVar(const std::string & var_name, unsigned int comp) const
281 : {
282 23314 : return getVarHelper<MooseVariableFieldBase>(var_name, comp);
283 : }
284 :
285 : MooseVariable *
286 2982 : Coupleable::getVar(const std::string & var_name, unsigned int comp)
287 : {
288 2982 : return const_cast<MooseVariable *>(getVarHelper<MooseVariable>(var_name, comp));
289 : }
290 :
291 : VectorMooseVariable *
292 557 : Coupleable::getVectorVar(const std::string & var_name, unsigned int comp)
293 : {
294 : auto * const var =
295 557 : const_cast<VectorMooseVariable *>(getVarHelper<VectorMooseVariable>(var_name, comp));
296 :
297 557 : if (_c_nodal && var && var->feType().family != LAGRANGE_VEC)
298 0 : mooseError(_c_name, ": Only LAGRANGE_VEC vector variables are defined at nodes");
299 :
300 557 : return var;
301 : }
302 :
303 : ArrayMooseVariable *
304 727 : Coupleable::getArrayVar(const std::string & var_name, unsigned int comp)
305 : {
306 727 : return const_cast<ArrayMooseVariable *>(getVarHelper<ArrayMooseVariable>(var_name, comp));
307 : }
308 :
309 : const MooseVariable *
310 25609 : Coupleable::getVar(const std::string & var_name, unsigned int comp) const
311 : {
312 25609 : return getVarHelper<MooseVariable>(var_name, comp);
313 : }
314 :
315 : const VectorMooseVariable *
316 1752 : Coupleable::getVectorVar(const std::string & var_name, unsigned int comp) const
317 : {
318 1752 : const auto * const var = getVarHelper<VectorMooseVariable>(var_name, comp);
319 :
320 1752 : if (_c_nodal && var && var->feType().family != LAGRANGE_VEC)
321 0 : mooseError(_c_name, ": Only LAGRANGE_VEC vector variables are defined at nodes");
322 :
323 1752 : return var;
324 : }
325 :
326 : const ArrayMooseVariable *
327 2259 : Coupleable::getArrayVar(const std::string & var_name, unsigned int comp) const
328 : {
329 2259 : return getVarHelper<ArrayMooseVariable>(var_name, comp);
330 : }
331 :
332 : const VariableValue *
333 1320 : Coupleable::getDefaultValue(const std::string & var_name, unsigned int comp) const
334 : {
335 : // make sure we don't access values that were not provided
336 1320 : checkComponent(_obj, comp, _c_parameters.numberDefaultCoupledValues(var_name), var_name);
337 :
338 1316 : auto default_value_it = _default_value.find(var_name);
339 1316 : if (default_value_it == _default_value.end())
340 : {
341 2568 : _default_value[var_name].emplace_back(std::make_unique<VariableValue>(
342 1284 : _coupleable_max_qps, _c_parameters.defaultCoupledValue(var_name, 0)));
343 1316 : for (unsigned int j = 1; j < _c_parameters.numberDefaultCoupledValues(var_name); ++j)
344 64 : _default_value[var_name].emplace_back(std::make_unique<VariableValue>(
345 64 : _coupleable_max_qps, _c_parameters.defaultCoupledValue(var_name, j)));
346 1284 : default_value_it = _default_value.find(var_name);
347 : }
348 :
349 1316 : const auto & default_value_vec = default_value_it->second;
350 1316 : const auto n_default_vals = default_value_vec.size();
351 1316 : if (comp >= n_default_vals)
352 0 : mooseError("Requested comp ",
353 : comp,
354 : " is equal to or greater than the number of default values ",
355 : n_default_vals);
356 2632 : return default_value_vec[comp].get();
357 : }
358 :
359 : const VectorVariableValue *
360 311 : Coupleable::getDefaultVectorValue(const std::string & var_name) const
361 : {
362 311 : auto default_value_it = _default_vector_value.find(var_name);
363 311 : if (default_value_it == _default_vector_value.end())
364 : {
365 311 : auto value = std::make_unique<VectorVariableValue>(_coupleable_max_qps, 0);
366 311 : bool already_warned = false;
367 311311 : for (unsigned int qp = 0; qp < _coupleable_max_qps; ++qp)
368 1244000 : for (const auto i : make_range(Moose::dim))
369 : {
370 : try
371 : {
372 933000 : (*value)[qp](i) = _c_parameters.defaultCoupledValue(var_name, i);
373 : }
374 11000 : catch (const std::out_of_range &)
375 : {
376 11000 : if (!already_warned)
377 11 : mooseWarning(
378 : "You supplied less than 3 arguments for the default vector value for variable ",
379 : var_name,
380 : ". Did you accidently leave something off? We are going to assign 0s, assuming "
381 : "this "
382 : "was intentional.");
383 11000 : already_warned = true;
384 11000 : (*value)[qp](i) = 0;
385 11000 : }
386 : }
387 311 : default_value_it =
388 311 : _default_vector_value.insert(std::make_pair(var_name, std::move(value))).first;
389 311 : }
390 :
391 622 : return default_value_it->second.get();
392 : }
393 :
394 : const ArrayVariableValue *
395 14 : Coupleable::getDefaultArrayValue(const std::string & var_name) const
396 : {
397 14 : auto default_value_it = _default_array_value.find(var_name);
398 14 : if (default_value_it == _default_array_value.end())
399 : {
400 14 : auto value = std::make_unique<ArrayVariableValue>(_coupleable_max_qps);
401 14014 : for (unsigned int qp = 0; qp < _coupleable_max_qps; ++qp)
402 : {
403 14000 : auto n = _c_parameters.numberDefaultCoupledValues(var_name);
404 14000 : (*value)[qp].resize(n);
405 42000 : for (unsigned int i = 0; i < n; ++i)
406 28000 : (*value)[qp](i) = _c_parameters.defaultCoupledValue(var_name, i);
407 : }
408 14 : default_value_it =
409 14 : _default_array_value.insert(std::make_pair(var_name, std::move(value))).first;
410 14 : }
411 :
412 28 : return default_value_it->second.get();
413 : }
414 :
415 : template <typename T>
416 : const T &
417 0 : Coupleable::getDefaultNodalValue(const std::string & var_name, unsigned int comp) const
418 : {
419 0 : auto && default_variable_value = getDefaultValue(var_name, comp);
420 0 : return *default_variable_value->data();
421 : }
422 :
423 : template <>
424 : const RealVectorValue &
425 0 : Coupleable::getDefaultNodalValue<RealVectorValue>(const std::string & var_name, unsigned int) const
426 : {
427 0 : auto && default_variable_value = getDefaultVectorValue(var_name);
428 0 : return *default_variable_value->data();
429 : }
430 :
431 : template <>
432 : const RealEigenVector &
433 0 : Coupleable::getDefaultNodalValue<RealEigenVector>(const std::string & var_name, unsigned int) const
434 : {
435 0 : auto && default_variable_value = getDefaultArrayValue(var_name);
436 0 : return *default_variable_value->data();
437 : }
438 :
439 : unsigned int
440 12278 : Coupleable::coupled(const std::string & var_name, unsigned int comp) const
441 : {
442 12278 : const auto * var = getFieldVar(var_name, comp);
443 12270 : if (!var)
444 : {
445 : mooseAssert(_optional_var_index.find(var_name) != _optional_var_index.end(),
446 : "optional var index for " << var_name << " does not exist!");
447 : // make sure we don't try to access default var ids that were not provided
448 142 : checkComponent(_obj, comp, _optional_var_index.at(var_name).size(), var_name);
449 138 : return _optional_var_index.at(var_name)[comp];
450 : }
451 12128 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
452 :
453 20770 : if (var->kind() == Moose::VAR_SOLVER &&
454 : // are we not an object that feeds into the nonlinear system?
455 8642 : (!_c_sys || _c_sys->varKind() != Moose::VAR_SOLVER ||
456 : // are we an object that impacts the nonlinear system and this variable is within our
457 : // nonlinear system?
458 4923 : var->sys().number() == _c_sys->number()))
459 8455 : return var->number();
460 : else
461 : // Avoid registering coupling to variables outside of our system (e.g. avoid potentially
462 : // creating bad Jacobians)
463 3673 : return std::numeric_limits<unsigned int>::max() - var->number();
464 : }
465 :
466 : template <>
467 : const GenericVariableValue<false> &
468 5854 : Coupleable::coupledGenericValue<false>(const std::string & var_name, unsigned int comp) const
469 : {
470 5854 : return coupledValue(var_name, comp);
471 : }
472 :
473 : template <>
474 : const GenericVariableValue<true> &
475 380 : Coupleable::coupledGenericValue<true>(const std::string & var_name, unsigned int comp) const
476 : {
477 380 : return adCoupledValue(var_name, comp);
478 : }
479 :
480 : template <>
481 : const GenericVectorVariableValue<false> &
482 261 : Coupleable::coupledGenericVectorValue<false>(const std::string & var_name, unsigned int comp) const
483 : {
484 261 : return coupledVectorValue(var_name, comp);
485 : }
486 :
487 : template <>
488 : const GenericVectorVariableValue<true> &
489 14 : Coupleable::coupledGenericVectorValue<true>(const std::string & var_name, unsigned int comp) const
490 : {
491 14 : return adCoupledVectorValue(var_name, comp);
492 : }
493 :
494 : const VariableValue &
495 57026 : Coupleable::coupledValue(const std::string & var_name, unsigned int comp) const
496 : {
497 57026 : const auto * const var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
498 57026 : if (!var)
499 1320 : return *getDefaultValue(var_name, comp);
500 55706 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
501 :
502 55706 : if (!_coupleable_neighbor)
503 : {
504 54198 : if (_c_nodal)
505 4353 : return (_c_is_implicit) ? var->dofValues() : var->dofValuesOld();
506 : else
507 49845 : return (_c_is_implicit) ? var->sln() : var->slnOld();
508 : }
509 : else
510 : {
511 1508 : if (_c_nodal)
512 0 : return (_c_is_implicit) ? var->dofValuesNeighbor() : var->dofValuesOldNeighbor();
513 : else
514 1508 : return (_c_is_implicit) ? var->slnNeighbor() : var->slnOldNeighbor();
515 : }
516 : }
517 :
518 : template <typename T>
519 : const typename OutputTools<T>::VariableValue &
520 1171 : Coupleable::vectorTagValueHelper(const std::string & var_names,
521 : const TagID tag,
522 : const unsigned int index) const
523 : {
524 1171 : const auto * const var = getVarHelper<MooseVariableField<T>>(var_names, index);
525 1171 : if (!var)
526 0 : mooseError(var_names, ": invalid variable name for coupledVectorTagValue");
527 1171 : checkFuncType(var_names, VarType::Ignore, FuncAge::Curr);
528 :
529 1171 : if (!_c_fe_problem.vectorTagExists(tag))
530 0 : mooseError("Attempting to couple to vector tag with ID ",
531 : tag,
532 : "in ",
533 0 : _c_name,
534 : ", but a vector tag with that ID does not exist");
535 :
536 1171 : const_cast<Coupleable *>(this)->addFEVariableCoupleableVectorTag(tag);
537 :
538 1171 : if (_c_nodal)
539 388 : return var->nodalVectorTagValue(tag);
540 : else
541 783 : return var->vectorTagValue(tag);
542 : }
543 :
544 : template <typename T>
545 : void
546 134 : Coupleable::requestStates(const std::string & var_name,
547 : const TagName & tag_name,
548 : const unsigned int comp)
549 : {
550 : auto var =
551 134 : const_cast<MooseVariableField<T> *>(getVarHelper<MooseVariableField<T>>(var_name, comp));
552 134 : if (!var)
553 0 : mooseError(var_name, ": invalid variable name for tag coupling");
554 :
555 134 : auto & var_sys = var->sys();
556 134 : if (tag_name == Moose::OLD_SOLUTION_TAG)
557 106 : var_sys.needSolutionState(1);
558 28 : else if (tag_name == Moose::OLDER_SOLUTION_TAG)
559 28 : var_sys.needSolutionState(2);
560 134 : }
561 :
562 : template <typename T>
563 : const typename OutputTools<T>::VariableValue &
564 1157 : Coupleable::vectorTagValueHelper(const std::string & var_names,
565 : const std::string & tag_param_name,
566 : const unsigned int index) const
567 : {
568 1157 : if (!_c_parameters.isParamValid(tag_param_name))
569 0 : mooseError("Tag name parameter '", tag_param_name, "' is invalid");
570 :
571 1157 : const TagName tag_name = MooseUtils::toUpper(_c_parameters.get<TagName>(tag_param_name));
572 :
573 1157 : const bool older_state_tag = _older_state_tags.count(tag_name);
574 1157 : if (older_state_tag)
575 : // We may need to add solution states and create vector tags
576 106 : const_cast<Coupleable *>(this)->requestStates<T>(var_names, tag_name, index);
577 :
578 1157 : if (!_c_fe_problem.vectorTagExists(tag_name))
579 0 : mooseError("Tagged vector with tag name '", tag_name, "' does not exist");
580 :
581 1157 : TagID tag = _c_fe_problem.getVectorTagID(tag_name);
582 2314 : return vectorTagValueHelper<T>(var_names, tag, index);
583 1157 : }
584 :
585 : template <>
586 : const GenericVariableValue<false> &
587 0 : Coupleable::coupledGenericDofValue<false>(const std::string & var_name, unsigned int comp) const
588 : {
589 0 : return coupledDofValues(var_name, comp);
590 : }
591 :
592 : template <>
593 : const GenericVariableValue<true> &
594 0 : Coupleable::coupledGenericDofValue<true>(const std::string & var_name, unsigned int comp) const
595 : {
596 0 : return adCoupledDofValues(var_name, comp);
597 : }
598 :
599 : const VariableValue &
600 14 : Coupleable::coupledValueLower(const std::string & var_name, const unsigned int comp) const
601 : {
602 14 : const auto * var = getVar(var_name, comp);
603 14 : if (!var)
604 0 : return *getDefaultValue(var_name, comp);
605 14 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
606 :
607 14 : if (_coupleable_neighbor)
608 0 : mooseError(_c_name, ":coupledValueLower cannot be called in a coupleable neighbor object");
609 :
610 14 : if (_c_nodal)
611 14 : return (_c_is_implicit) ? var->dofValues() : var->dofValuesOld();
612 : else
613 0 : return (_c_is_implicit) ? var->slnLower() : var->slnLowerOld();
614 : }
615 :
616 : const VariableValue &
617 14 : Coupleable::coupledVectorTagValue(const std::string & var_names,
618 : TagID tag,
619 : unsigned int index) const
620 : {
621 14 : return vectorTagValueHelper<Real>(var_names, tag, index);
622 : }
623 :
624 : const VariableValue &
625 1107 : Coupleable::coupledVectorTagValue(const std::string & var_names,
626 : const std::string & tag_name,
627 : unsigned int index) const
628 : {
629 1107 : return vectorTagValueHelper<Real>(var_names, tag_name, index);
630 : }
631 :
632 : const ArrayVariableValue &
633 0 : Coupleable::coupledVectorTagArrayValue(const std::string & var_names,
634 : TagID tag,
635 : unsigned int index) const
636 : {
637 0 : return vectorTagValueHelper<RealEigenVector>(var_names, tag, index);
638 : }
639 :
640 : const ArrayVariableValue &
641 50 : Coupleable::coupledVectorTagArrayValue(const std::string & var_names,
642 : const std::string & tag_name,
643 : unsigned int index) const
644 : {
645 50 : return vectorTagValueHelper<RealEigenVector>(var_names, tag_name, index);
646 : }
647 :
648 : const VariableGradient &
649 0 : Coupleable::coupledVectorTagGradient(const std::string & var_names,
650 : TagID tag,
651 : unsigned int index) const
652 : {
653 0 : const auto * var = getVar(var_names, index);
654 0 : if (!var)
655 0 : mooseError(var_names, ": invalid variable name for coupledVectorTagGradient");
656 0 : checkFuncType(var_names, VarType::Ignore, FuncAge::Curr);
657 :
658 0 : if (!_c_fe_problem.vectorTagExists(tag))
659 0 : mooseError("Attempting to couple to vector tag with ID ",
660 : tag,
661 : "in ",
662 0 : _c_name,
663 : ", but a vector tag with that ID does not exist");
664 :
665 0 : const_cast<Coupleable *>(this)->addFEVariableCoupleableVectorTag(tag);
666 :
667 0 : return var->vectorTagGradient(tag);
668 : }
669 :
670 : const VariableGradient &
671 0 : Coupleable::coupledVectorTagGradient(const std::string & var_names,
672 : const std::string & tag_name,
673 : unsigned int index) const
674 : {
675 0 : if (!_c_parameters.isParamValid(tag_name))
676 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
677 :
678 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
679 0 : if (!_c_fe_problem.vectorTagExists(tagname))
680 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
681 :
682 0 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
683 0 : return coupledVectorTagGradient(var_names, tag, index);
684 0 : }
685 :
686 : const ArrayVariableGradient &
687 28 : Coupleable::coupledVectorTagArrayGradient(const std::string & var_names,
688 : TagID tag,
689 : unsigned int index) const
690 : {
691 28 : const auto * var = getArrayVar(var_names, index);
692 28 : if (!var)
693 0 : mooseError(var_names, ": invalid variable name for coupledVectorTagArrayGradient");
694 28 : checkFuncType(var_names, VarType::Ignore, FuncAge::Curr);
695 :
696 28 : if (!_c_fe_problem.vectorTagExists(tag))
697 0 : mooseError("Attempting to couple to vector tag with ID ",
698 : tag,
699 : "in ",
700 0 : _c_name,
701 : ", but a vector tag with that ID does not exist");
702 :
703 28 : const_cast<Coupleable *>(this)->addFEVariableCoupleableVectorTag(tag);
704 :
705 28 : return var->vectorTagGradient(tag);
706 : }
707 :
708 : const ArrayVariableGradient &
709 28 : Coupleable::coupledVectorTagArrayGradient(const std::string & var_names,
710 : const std::string & tag_name,
711 : unsigned int index) const
712 : {
713 28 : if (!_c_parameters.isParamValid(tag_name))
714 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
715 :
716 28 : TagName tagname = _c_parameters.get<TagName>(tag_name);
717 28 : if (!_c_fe_problem.vectorTagExists(tagname))
718 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
719 :
720 28 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
721 56 : return coupledVectorTagArrayGradient(var_names, tag, index);
722 28 : }
723 :
724 : template <typename T>
725 : const typename OutputTools<T>::VariableValue &
726 148 : Coupleable::vectorTagDofValueHelper(const std::string & var_name,
727 : const TagID tag,
728 : const unsigned int comp) const
729 : {
730 148 : const auto * var = getVarHelper<MooseVariableField<T>>(var_name, comp);
731 148 : if (!var)
732 0 : mooseError(var_name, ": invalid variable name for coupledVectorTagDofValue");
733 148 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
734 :
735 148 : const_cast<Coupleable *>(this)->addFEVariableCoupleableVectorTag(tag);
736 :
737 148 : return var->vectorTagDofValue(tag);
738 : }
739 :
740 : template <typename T>
741 : const typename OutputTools<T>::VariableValue &
742 148 : Coupleable::vectorTagDofValueHelper(const std::string & var_name,
743 : const std::string & tag_param_name,
744 : const unsigned int comp) const
745 : {
746 148 : if (!_c_parameters.isParamValid(tag_param_name))
747 0 : mooseError("Tag name parameter '", tag_param_name, "' is invalid");
748 :
749 148 : const TagName tag_name = MooseUtils::toUpper(_c_parameters.get<TagName>(tag_param_name));
750 :
751 148 : const bool older_state_tag = _older_state_tags.count(tag_name);
752 148 : if (older_state_tag)
753 : // We may need to add solution states and create vector tags
754 28 : const_cast<Coupleable *>(this)->requestStates<T>(var_name, tag_name, comp);
755 :
756 148 : if (!_c_fe_problem.vectorTagExists(tag_name))
757 0 : mooseError("Tagged vector with tag name '", tag_name, "' does not exist");
758 :
759 148 : TagID tag = _c_fe_problem.getVectorTagID(tag_name);
760 :
761 296 : return vectorTagDofValueHelper<T>(var_name, tag, comp);
762 148 : }
763 :
764 : const VariableValue &
765 0 : Coupleable::coupledVectorTagDofValue(const std::string & var_name,
766 : TagID tag,
767 : unsigned int comp) const
768 : {
769 0 : return vectorTagDofValueHelper<Real>(var_name, tag, comp);
770 : }
771 :
772 : const VariableValue &
773 134 : Coupleable::coupledVectorTagDofValue(const std::string & var_name,
774 : const std::string & tag_name,
775 : unsigned int comp) const
776 : {
777 134 : return vectorTagDofValueHelper<Real>(var_name, tag_name, comp);
778 : }
779 :
780 : const ArrayVariableValue &
781 14 : Coupleable::coupledVectorTagArrayDofValue(const std::string & var_name,
782 : const std::string & tag_name,
783 : unsigned int comp) const
784 : {
785 14 : return vectorTagDofValueHelper<RealEigenVector>(var_name, tag_name, comp);
786 : }
787 :
788 : const VariableValue &
789 194 : Coupleable::coupledMatrixTagValue(const std::string & var_names,
790 : TagID tag,
791 : unsigned int index) const
792 : {
793 194 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_names, index);
794 194 : if (!var)
795 0 : mooseError(var_names, ": invalid variable name for coupledMatrixTagValue");
796 194 : checkFuncType(var_names, VarType::Ignore, FuncAge::Curr);
797 :
798 194 : const_cast<Coupleable *>(this)->addFEVariableCoupleableMatrixTag(tag);
799 :
800 194 : if (_c_nodal)
801 140 : return var->nodalMatrixTagValue(tag);
802 54 : return var->matrixTagValue(tag);
803 : }
804 :
805 : const VariableValue &
806 194 : Coupleable::coupledMatrixTagValue(const std::string & var_names,
807 : const std::string & tag_name,
808 : unsigned int index) const
809 : {
810 194 : if (!_c_parameters.isParamValid(tag_name))
811 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
812 :
813 194 : TagName tagname = _c_parameters.get<TagName>(tag_name);
814 194 : if (!_c_fe_problem.matrixTagExists(tagname))
815 0 : mooseError("Matrix tag name '", tagname, "' does not exist");
816 :
817 194 : TagID tag = _c_fe_problem.getMatrixTagID(tagname);
818 388 : return coupledMatrixTagValue(var_names, tag, index);
819 194 : }
820 :
821 : const VectorVariableValue &
822 935 : Coupleable::coupledVectorValue(const std::string & var_name, unsigned int comp) const
823 : {
824 935 : const auto * var = getVectorVar(var_name, comp);
825 935 : if (!var)
826 311 : return *getDefaultVectorValue(var_name);
827 624 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
828 :
829 624 : if (!_coupleable_neighbor)
830 : {
831 596 : if (_c_nodal)
832 70 : return _c_is_implicit ? var->nodalValueArray() : var->nodalValueOldArray();
833 : else
834 526 : return _c_is_implicit ? var->sln() : var->slnOld();
835 : }
836 : else
837 : {
838 28 : if (_c_nodal)
839 : // Since this is at a node, I don't feel like there should be any "neighbor" logic
840 0 : return _c_is_implicit ? var->nodalValueArray() : var->nodalValueOldArray();
841 : else
842 28 : return _c_is_implicit ? var->slnNeighbor() : var->slnOldNeighbor();
843 : }
844 : }
845 :
846 : const ArrayVariableValue &
847 1217 : Coupleable::coupledArrayValue(const std::string & var_name, unsigned int comp) const
848 : {
849 1217 : const auto * var = getArrayVar(var_name, comp);
850 1217 : if (!var)
851 14 : return *getDefaultArrayValue(var_name);
852 1203 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
853 :
854 1203 : if (!_coupleable_neighbor)
855 : {
856 1203 : if (_c_nodal)
857 92 : return (_c_is_implicit) ? var->dofValues() : var->dofValuesOld();
858 1111 : return (_c_is_implicit) ? var->sln() : var->slnOld();
859 : }
860 : else
861 : {
862 0 : if (_c_nodal)
863 0 : return (_c_is_implicit) ? var->dofValuesNeighbor() : var->dofValuesOldNeighbor();
864 0 : return (_c_is_implicit) ? var->slnNeighbor() : var->slnOldNeighbor();
865 : }
866 : }
867 :
868 : std::vector<const ArrayVariableValue *>
869 18 : Coupleable::coupledArrayValues(const std::string & var_name) const
870 : {
871 36 : auto func = [this, &var_name](unsigned int comp) { return &coupledArrayValue(var_name, comp); };
872 36 : return coupledVectorHelper<const ArrayVariableValue *>(var_name, func);
873 : }
874 :
875 : MooseWritableVariable &
876 282 : Coupleable::writableVariable(const std::string & var_name, unsigned int comp)
877 : {
878 282 : auto * var = getVarHelper<MooseWritableVariable>(var_name, comp);
879 :
880 278 : const auto * aux = dynamic_cast<const AuxKernel *>(this);
881 278 : const auto * euo = dynamic_cast<const ElementUserObject *>(this);
882 278 : const auto * nuo = dynamic_cast<const NodalUserObject *>(this);
883 278 : const auto * nfc = dynamic_cast<const NodeFaceConstraint *>(this);
884 278 : const auto * nec = dynamic_cast<const NodeElemConstraintBase *>(this);
885 :
886 278 : if (!aux && !euo && !nuo && !nfc && !nec)
887 4 : mooseError("writableVariable() can only be called from AuxKernels, ElementUserObjects, "
888 : "NodalUserObjects, NodeFaceConstraints, or NodeElemConstraints. '",
889 4 : _obj->name(),
890 : "' is none of those.");
891 :
892 274 : if (aux && !aux->isNodal() && var->isNodal())
893 4 : mooseError("The elemental AuxKernel '",
894 4 : _obj->name(),
895 : "' cannot obtain a writable reference to the nodal variable '",
896 4 : var->name(),
897 : "'.");
898 270 : if (euo && var->isNodal())
899 4 : mooseError("The ElementUserObject '",
900 4 : _obj->name(),
901 : "' cannot obtain a writable reference to the nodal variable '",
902 4 : var->name(),
903 : "'.");
904 :
905 : // make sure only one object can access a variable
906 266 : checkWritableVar(var);
907 :
908 254 : return *var;
909 : }
910 :
911 : VariableValue &
912 44 : Coupleable::writableCoupledValue(const std::string & var_name, unsigned int comp)
913 : {
914 44 : mooseDeprecated("Coupleable::writableCoupledValue is deprecated, please use "
915 : "Coupleable::writableVariable instead. ");
916 :
917 : // check if the variable exists
918 44 : auto * const var = getVar(var_name, comp);
919 44 : if (!var)
920 4 : mooseError(
921 : "Unable to create a writable reference for '", var_name, "', is it a constant expression?");
922 :
923 : // is the requested variable an AuxiliaryVariable?
924 40 : if (!_c_fe_problem.getAuxiliarySystem().hasVariable(var->name()))
925 4 : mooseError(
926 4 : "'", var->name(), "' must be an auxiliary variable in Coupleable::writableCoupledValue");
927 :
928 : // check that the variable type (elemental/nodal) is compatible with the object type
929 36 : const auto * aux = dynamic_cast<const AuxKernel *>(this);
930 :
931 36 : if (!aux)
932 4 : mooseError("writableCoupledValue() can only be called from AuxKernels, but '",
933 4 : _obj->name(),
934 : "' is not an AuxKernel.");
935 :
936 32 : if (!aux->isNodal() && var->isNodal())
937 4 : mooseError("The elemental AuxKernel '",
938 4 : _obj->name(),
939 : "' cannot obtain a writable reference to the nodal variable '",
940 4 : var->name(),
941 : "'.");
942 :
943 : // make sure only one object can access a variable
944 28 : checkWritableVar(var);
945 :
946 28 : return const_cast<VariableValue &>(coupledValue(var_name, comp));
947 : }
948 :
949 : void
950 294 : Coupleable::checkWritableVar(MooseWritableVariable * var)
951 : {
952 : // check domain restrictions for compatibility
953 294 : const auto * br = dynamic_cast<const BlockRestrictable *>(this);
954 294 : const auto * nfc = dynamic_cast<const NodeFaceConstraint *>(this);
955 :
956 294 : if (br && !var->hasBlocks(br->blockIDs()))
957 4 : mooseError("The variable '",
958 4 : var->name(),
959 : "' must be defined on all blocks '",
960 4 : _obj->name(),
961 : "' is defined on.");
962 :
963 290 : if (nfc && !var->hasBlocks(nfc->getSecondaryConnectedBlocks()))
964 0 : mooseError("The variable '",
965 0 : var->name(),
966 : " must be defined on all blocks '",
967 0 : _obj->name(),
968 : "'s secondary surface is defined on.");
969 :
970 : // make sure only one object can access a variable
971 1080 : for (const auto & ci : _obj->getMooseApp().getInterfaceObjects<Coupleable>())
972 798 : if (ci != this && ci->_writable_coupled_variables[_c_tid].count(var))
973 : {
974 : // if both this and ci are block restrictable then we check if the block restrictions
975 : // are not overlapping. If they don't we permit the call.
976 41 : const auto * br_other = dynamic_cast<const BlockRestrictable *>(ci);
977 78 : if (br && br_other && br->blockRestricted() && br_other->blockRestricted() &&
978 37 : !MooseUtils::setsIntersect(br->blockIDs(), br_other->blockIDs()))
979 33 : continue;
980 8 : else if (nfc)
981 0 : continue;
982 :
983 8 : mooseError("'",
984 8 : ci->_obj->name(),
985 : "' already obtained a writable reference to '",
986 8 : var->name(),
987 : "'. Only one object can obtain such a reference per variable and subdomain in a "
988 : "simulation.");
989 : }
990 :
991 : // var is unique across threads, so we could forego having a separate set per thread, but we
992 : // need quick access to the list of all variables that need to be inserted into the solution
993 : // vector by a given thread.
994 :
995 282 : _writable_coupled_variables[_c_tid].insert(var);
996 282 : }
997 :
998 : const VariableValue &
999 12250 : Coupleable::coupledValueOld(const std::string & var_name, unsigned int comp) const
1000 : {
1001 12250 : const auto * var = getVar(var_name, comp);
1002 12250 : if (!var)
1003 0 : return *getDefaultValue(var_name, comp);
1004 12250 : checkFuncType(var_name, VarType::Ignore, FuncAge::Old);
1005 :
1006 12250 : if (!_coupleable_neighbor)
1007 : {
1008 8279 : if (_c_nodal)
1009 96 : return (_c_is_implicit) ? var->dofValuesOld() : var->dofValuesOlder();
1010 8183 : return (_c_is_implicit) ? var->slnOld() : var->slnOlder();
1011 : }
1012 : else
1013 : {
1014 3971 : if (_c_nodal)
1015 0 : return (_c_is_implicit) ? var->dofValuesOldNeighbor() : var->dofValuesOlderNeighbor();
1016 3971 : return (_c_is_implicit) ? var->slnOldNeighbor() : var->slnOlderNeighbor();
1017 : }
1018 : }
1019 :
1020 : const VariableValue &
1021 11928 : Coupleable::coupledValueOlder(const std::string & var_name, unsigned int comp) const
1022 : {
1023 11928 : const auto * var = getVar(var_name, comp);
1024 11928 : if (!var)
1025 0 : return *getDefaultValue(var_name, comp);
1026 11928 : checkFuncType(var_name, VarType::Ignore, FuncAge::Older);
1027 :
1028 11928 : if (!_coupleable_neighbor)
1029 : {
1030 7980 : if (_c_nodal)
1031 42 : return var->dofValuesOlder();
1032 7938 : return var->slnOlder();
1033 : }
1034 : else
1035 : {
1036 3948 : if (_c_nodal)
1037 0 : return var->dofValuesOlderNeighbor();
1038 3948 : return var->slnOlderNeighbor();
1039 : }
1040 : }
1041 :
1042 : const VariableValue &
1043 82 : Coupleable::coupledValuePreviousNL(const std::string & var_name, unsigned int comp) const
1044 : {
1045 82 : const auto * var = getVar(var_name, comp);
1046 82 : if (!var)
1047 0 : return *getDefaultValue(var_name, comp);
1048 82 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
1049 :
1050 82 : _c_fe_problem.needsPreviousNewtonIteration(true);
1051 82 : if (!_coupleable_neighbor)
1052 : {
1053 82 : if (_c_nodal)
1054 0 : return var->dofValuesPreviousNL();
1055 82 : return var->slnPreviousNL();
1056 : }
1057 : else
1058 : {
1059 0 : if (_c_nodal)
1060 0 : return var->dofValuesPreviousNLNeighbor();
1061 0 : return var->slnPreviousNLNeighbor();
1062 : }
1063 : }
1064 :
1065 : const VectorVariableValue &
1066 56 : Coupleable::coupledVectorValueOld(const std::string & var_name, unsigned int comp) const
1067 : {
1068 56 : const auto * var = getVectorVar(var_name, comp);
1069 56 : if (!var)
1070 0 : return *getDefaultVectorValue(var_name);
1071 56 : checkFuncType(var_name, VarType::Ignore, FuncAge::Old);
1072 :
1073 56 : if (_c_nodal)
1074 28 : return (_c_is_implicit) ? var->nodalValueOldArray() : var->nodalValueOlderArray();
1075 28 : if (!_coupleable_neighbor)
1076 28 : return (_c_is_implicit) ? var->slnOld() : var->slnOlder();
1077 0 : return (_c_is_implicit) ? var->slnOldNeighbor() : var->slnOlderNeighbor();
1078 : }
1079 :
1080 : const VectorVariableValue &
1081 0 : Coupleable::coupledVectorValueOlder(const std::string & var_name, unsigned int comp) const
1082 : {
1083 0 : const auto * var = getVectorVar(var_name, comp);
1084 0 : if (!var)
1085 0 : return *getDefaultVectorValue(var_name);
1086 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Older);
1087 :
1088 0 : if (!_coupleable_neighbor)
1089 0 : return var->slnOlder();
1090 0 : return var->slnOlderNeighbor();
1091 : }
1092 :
1093 : const ArrayVariableValue &
1094 0 : Coupleable::coupledArrayValueOld(const std::string & var_name, unsigned int comp) const
1095 : {
1096 0 : const auto * var = getArrayVar(var_name, comp);
1097 0 : if (!var)
1098 0 : return *getDefaultArrayValue(var_name);
1099 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Old);
1100 :
1101 0 : if (!_coupleable_neighbor)
1102 : {
1103 0 : if (_c_nodal)
1104 0 : return (_c_is_implicit) ? var->dofValuesOld() : var->dofValuesOlder();
1105 0 : return (_c_is_implicit) ? var->slnOld() : var->slnOlder();
1106 : }
1107 : else
1108 : {
1109 0 : if (_c_nodal)
1110 0 : return (_c_is_implicit) ? var->dofValuesOldNeighbor() : var->dofValuesOlderNeighbor();
1111 0 : return (_c_is_implicit) ? var->slnOldNeighbor() : var->slnOlderNeighbor();
1112 : }
1113 : }
1114 :
1115 : const ArrayVariableValue &
1116 0 : Coupleable::coupledArrayValueOlder(const std::string & var_name, unsigned int comp) const
1117 : {
1118 0 : const auto * var = getArrayVar(var_name, comp);
1119 0 : if (!var)
1120 0 : return *getDefaultArrayValue(var_name);
1121 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Older);
1122 :
1123 0 : if (!_coupleable_neighbor)
1124 : {
1125 0 : if (_c_nodal)
1126 0 : return var->dofValuesOlder();
1127 0 : return var->slnOlder();
1128 : }
1129 : else
1130 : {
1131 0 : if (_c_nodal)
1132 0 : return var->dofValuesOlderNeighbor();
1133 0 : return var->slnOlderNeighbor();
1134 : }
1135 : }
1136 :
1137 : const VariableValue &
1138 296 : Coupleable::coupledDot(const std::string & var_name, unsigned int comp) const
1139 : {
1140 296 : const auto * var = getVar(var_name, comp);
1141 296 : if (!var)
1142 : {
1143 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1144 0 : return _default_value_zero;
1145 : }
1146 296 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1147 :
1148 292 : if (!_coupleable_neighbor)
1149 : {
1150 264 : if (_c_nodal)
1151 42 : return var->dofValuesDot();
1152 222 : return var->uDot();
1153 : }
1154 : else
1155 : {
1156 28 : if (_c_nodal)
1157 0 : return var->dofValuesDotNeighbor();
1158 28 : return var->uDotNeighbor();
1159 : }
1160 : }
1161 :
1162 : const VariableValue &
1163 181 : Coupleable::coupledDotDot(const std::string & var_name, unsigned int comp) const
1164 : {
1165 181 : const auto * var = getVar(var_name, comp);
1166 181 : if (!var)
1167 : {
1168 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1169 0 : return _default_value_zero;
1170 : }
1171 181 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1172 :
1173 181 : if (!_coupleable_neighbor)
1174 : {
1175 153 : if (_c_nodal)
1176 4 : return var->dofValuesDotDot();
1177 149 : return var->uDotDot();
1178 : }
1179 : else
1180 : {
1181 28 : if (_c_nodal)
1182 0 : return var->dofValuesDotDotNeighbor();
1183 28 : return var->uDotDotNeighbor();
1184 : }
1185 : }
1186 :
1187 : template <>
1188 : const GenericVariableValue<false> &
1189 84 : Coupleable::coupledGenericDotDot<false>(const std::string & var_name, unsigned int comp) const
1190 : {
1191 84 : return coupledDotDot(var_name, comp);
1192 : }
1193 :
1194 : template <>
1195 : const GenericVariableValue<true> &
1196 42 : Coupleable::coupledGenericDotDot<true>(const std::string & var_name, unsigned int comp) const
1197 : {
1198 42 : return adCoupledDotDot(var_name, comp);
1199 : }
1200 :
1201 : const VariableValue &
1202 0 : Coupleable::coupledDotOld(const std::string & var_name, unsigned int comp) const
1203 : {
1204 0 : const auto * var = getVar(var_name, comp);
1205 0 : if (!var)
1206 : {
1207 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1208 0 : return _default_value_zero;
1209 : }
1210 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1211 :
1212 0 : if (!_coupleable_neighbor)
1213 : {
1214 0 : if (_c_nodal)
1215 0 : return var->dofValuesDotOld();
1216 0 : return var->uDotOld();
1217 : }
1218 : else
1219 : {
1220 0 : if (_c_nodal)
1221 0 : return var->dofValuesDotOldNeighbor();
1222 0 : return var->uDotOldNeighbor();
1223 : }
1224 : }
1225 :
1226 : const VariableValue &
1227 0 : Coupleable::coupledDotDotOld(const std::string & var_name, unsigned int comp) const
1228 : {
1229 0 : const auto * var = getVar(var_name, comp);
1230 0 : if (!var)
1231 : {
1232 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1233 0 : return _default_value_zero;
1234 : }
1235 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1236 :
1237 0 : if (!_coupleable_neighbor)
1238 : {
1239 0 : if (_c_nodal)
1240 0 : return var->dofValuesDotDotOld();
1241 0 : return var->uDotDotOld();
1242 : }
1243 : else
1244 : {
1245 0 : if (_c_nodal)
1246 0 : return var->dofValuesDotDotOldNeighbor();
1247 0 : return var->uDotDotOldNeighbor();
1248 : }
1249 : }
1250 :
1251 : const VectorVariableValue &
1252 123 : Coupleable::coupledVectorDot(const std::string & var_name, unsigned int comp) const
1253 : {
1254 123 : const auto * var = getVectorVar(var_name, comp);
1255 123 : if (!var)
1256 : {
1257 0 : _default_vector_value_zero.resize(_coupleable_max_qps);
1258 0 : return _default_vector_value_zero;
1259 : }
1260 123 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1261 :
1262 123 : if (!_coupleable_neighbor)
1263 95 : return var->uDot();
1264 28 : return var->uDotNeighbor();
1265 : }
1266 :
1267 : const VectorVariableValue &
1268 42 : Coupleable::coupledVectorDotDot(const std::string & var_name, unsigned int comp) const
1269 : {
1270 42 : const auto * var = getVectorVar(var_name, comp);
1271 42 : if (!var)
1272 : {
1273 0 : _default_vector_value_zero.resize(_coupleable_max_qps);
1274 0 : return _default_vector_value_zero;
1275 : }
1276 42 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1277 :
1278 42 : if (!_coupleable_neighbor)
1279 28 : return var->uDotDot();
1280 14 : return var->uDotDotNeighbor();
1281 : }
1282 :
1283 : const VectorVariableValue &
1284 0 : Coupleable::coupledVectorDotOld(const std::string & var_name, unsigned int comp) const
1285 : {
1286 0 : const auto * var = getVectorVar(var_name, comp);
1287 0 : if (!var)
1288 : {
1289 0 : _default_vector_value_zero.resize(_coupleable_max_qps);
1290 0 : return _default_vector_value_zero;
1291 : }
1292 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1293 :
1294 0 : if (!_coupleable_neighbor)
1295 0 : return var->uDotOld();
1296 0 : return var->uDotOldNeighbor();
1297 : }
1298 :
1299 : const VectorVariableValue &
1300 0 : Coupleable::coupledVectorDotDotOld(const std::string & var_name, unsigned int comp) const
1301 : {
1302 0 : const auto * var = getVectorVar(var_name, comp);
1303 0 : if (!var)
1304 : {
1305 0 : _default_vector_value_zero.resize(_coupleable_max_qps);
1306 0 : return _default_vector_value_zero;
1307 : }
1308 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1309 :
1310 0 : if (!_coupleable_neighbor)
1311 0 : return var->uDotDotOld();
1312 0 : return var->uDotDotOldNeighbor();
1313 : }
1314 :
1315 : const VariableValue &
1316 123 : Coupleable::coupledVectorDotDu(const std::string & var_name, unsigned int comp) const
1317 : {
1318 123 : const auto * var = getVectorVar(var_name, comp);
1319 123 : if (!var)
1320 : {
1321 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1322 0 : return _default_value_zero;
1323 : }
1324 123 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1325 :
1326 123 : if (!_coupleable_neighbor)
1327 95 : return var->duDotDu();
1328 28 : return var->duDotDuNeighbor();
1329 : }
1330 :
1331 : const VariableValue &
1332 42 : Coupleable::coupledVectorDotDotDu(const std::string & var_name, unsigned int comp) const
1333 : {
1334 42 : const auto * var = getVectorVar(var_name, comp);
1335 42 : if (!var)
1336 : {
1337 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1338 0 : return _default_value_zero;
1339 : }
1340 42 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1341 :
1342 42 : if (!_coupleable_neighbor)
1343 28 : return var->duDotDotDu();
1344 14 : return var->duDotDotDuNeighbor();
1345 : }
1346 :
1347 : const ArrayVariableValue &
1348 47 : Coupleable::coupledArrayDot(const std::string & var_name, unsigned int comp) const
1349 : {
1350 47 : const auto * var = getArrayVar(var_name, comp);
1351 47 : if (!var)
1352 0 : return _default_array_value_zero;
1353 47 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1354 :
1355 47 : if (!_coupleable_neighbor)
1356 : {
1357 47 : if (_c_nodal)
1358 0 : return var->dofValuesDot();
1359 47 : return var->uDot();
1360 : }
1361 : else
1362 : {
1363 0 : if (_c_nodal)
1364 0 : return var->dofValuesDotNeighbor();
1365 0 : return var->uDotNeighbor();
1366 : }
1367 : }
1368 :
1369 : const ArrayVariableValue &
1370 0 : Coupleable::coupledArrayDotDot(const std::string & var_name, unsigned int comp) const
1371 : {
1372 0 : const auto * var = getArrayVar(var_name, comp);
1373 0 : if (!var)
1374 0 : return _default_array_value_zero;
1375 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1376 :
1377 0 : if (!_coupleable_neighbor)
1378 : {
1379 0 : if (_c_nodal)
1380 0 : return var->dofValuesDotDot();
1381 0 : return var->uDotDot();
1382 : }
1383 : else
1384 : {
1385 0 : if (_c_nodal)
1386 0 : return var->dofValuesDotDotNeighbor();
1387 0 : return var->uDotDotNeighbor();
1388 : }
1389 : }
1390 :
1391 : const ArrayVariableValue &
1392 0 : Coupleable::coupledArrayDotOld(const std::string & var_name, unsigned int comp) const
1393 : {
1394 0 : const auto * var = getArrayVar(var_name, comp);
1395 0 : if (!var)
1396 0 : return _default_array_value_zero;
1397 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1398 :
1399 0 : if (!_coupleable_neighbor)
1400 : {
1401 0 : if (_c_nodal)
1402 0 : return var->dofValuesDotOld();
1403 0 : return var->uDotOld();
1404 : }
1405 : else
1406 : {
1407 0 : if (_c_nodal)
1408 0 : return var->dofValuesDotOldNeighbor();
1409 0 : return var->uDotOldNeighbor();
1410 : }
1411 : }
1412 :
1413 : const ArrayVariableValue &
1414 0 : Coupleable::coupledArrayDotDotOld(const std::string & var_name, unsigned int comp) const
1415 : {
1416 0 : const auto * var = getArrayVar(var_name, comp);
1417 0 : if (!var)
1418 0 : return _default_array_value_zero;
1419 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
1420 :
1421 0 : if (!_coupleable_neighbor)
1422 : {
1423 0 : if (_c_nodal)
1424 0 : return var->dofValuesDotDotOld();
1425 0 : return var->uDotDotOld();
1426 : }
1427 : else
1428 : {
1429 0 : if (_c_nodal)
1430 0 : return var->dofValuesDotDotOldNeighbor();
1431 0 : return var->uDotDotOldNeighbor();
1432 : }
1433 : }
1434 :
1435 : const VariableValue &
1436 179 : Coupleable::coupledDotDu(const std::string & var_name, unsigned int comp) const
1437 : {
1438 179 : const auto * var = getVar(var_name, comp);
1439 179 : if (!var)
1440 : {
1441 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1442 0 : return _default_value_zero;
1443 : }
1444 179 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1445 :
1446 179 : if (!_coupleable_neighbor)
1447 : {
1448 137 : if (_c_nodal)
1449 0 : return var->dofValuesDuDotDu();
1450 137 : return var->duDotDu();
1451 : }
1452 : else
1453 : {
1454 42 : if (_c_nodal)
1455 0 : return var->dofValuesDuDotDuNeighbor();
1456 42 : return var->duDotDuNeighbor();
1457 : }
1458 : }
1459 :
1460 : const VariableValue &
1461 126 : Coupleable::coupledDotDotDu(const std::string & var_name, unsigned int comp) const
1462 : {
1463 126 : const auto * var = getVar(var_name, comp);
1464 126 : if (!var)
1465 : {
1466 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1467 0 : return _default_value_zero;
1468 : }
1469 126 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1470 :
1471 126 : if (!_coupleable_neighbor)
1472 : {
1473 84 : if (_c_nodal)
1474 0 : return var->dofValuesDuDotDotDu();
1475 84 : return var->duDotDotDu();
1476 : }
1477 : else
1478 : {
1479 42 : if (_c_nodal)
1480 0 : return var->dofValuesDuDotDotDuNeighbor();
1481 42 : return var->duDotDotDuNeighbor();
1482 : }
1483 : }
1484 :
1485 : const VariableValue &
1486 47 : Coupleable::coupledArrayDotDu(const std::string & var_name, unsigned int comp) const
1487 : {
1488 47 : const auto * const var = getArrayVar(var_name, comp);
1489 47 : if (!var)
1490 : {
1491 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1492 0 : return _default_value_zero;
1493 : }
1494 47 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1495 :
1496 47 : if (!_coupleable_neighbor)
1497 : {
1498 47 : if (_c_nodal)
1499 0 : return var->dofValuesDuDotDu();
1500 47 : return var->duDotDu();
1501 : }
1502 : else
1503 : {
1504 0 : if (_c_nodal)
1505 0 : return var->dofValuesDuDotDuNeighbor();
1506 0 : return var->duDotDuNeighbor();
1507 : }
1508 : }
1509 :
1510 : const VariableGradient &
1511 32480 : Coupleable::coupledGradient(const std::string & var_name, unsigned int comp) const
1512 : {
1513 32480 : const auto * const var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
1514 32476 : if (!var)
1515 : {
1516 38 : _default_gradient.resize(_coupleable_max_qps);
1517 38 : return _default_gradient;
1518 : }
1519 32438 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1520 :
1521 32438 : if (!_coupleable_neighbor)
1522 32424 : return (_c_is_implicit) ? var->gradSln() : var->gradSlnOld();
1523 14 : return (_c_is_implicit) ? var->gradSlnNeighbor() : var->gradSlnOldNeighbor();
1524 : }
1525 :
1526 : const VariableGradient &
1527 4 : Coupleable::coupledGradientOld(const std::string & var_name, unsigned int comp) const
1528 : {
1529 4 : const auto * var = getVar(var_name, comp);
1530 4 : if (!var)
1531 : {
1532 0 : _default_gradient.resize(_coupleable_max_qps);
1533 0 : return _default_gradient;
1534 : }
1535 4 : checkFuncType(var_name, VarType::Gradient, FuncAge::Old);
1536 :
1537 0 : if (!_coupleable_neighbor)
1538 0 : return (_c_is_implicit) ? var->gradSlnOld() : var->gradSlnOlder();
1539 0 : return (_c_is_implicit) ? var->gradSlnOldNeighbor() : var->gradSlnOlderNeighbor();
1540 : }
1541 :
1542 : const VariableGradient &
1543 0 : Coupleable::coupledGradientOlder(const std::string & var_name, unsigned int comp) const
1544 : {
1545 0 : const auto * var = getVar(var_name, comp);
1546 0 : if (!var)
1547 : {
1548 0 : _default_gradient.resize(_coupleable_max_qps);
1549 0 : return _default_gradient;
1550 : }
1551 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1552 :
1553 0 : if (!_coupleable_neighbor)
1554 0 : return var->gradSlnOlder();
1555 0 : return var->gradSlnOlderNeighbor();
1556 : }
1557 :
1558 : const VariableGradient &
1559 0 : Coupleable::coupledGradientPreviousNL(const std::string & var_name, unsigned int comp) const
1560 : {
1561 0 : const auto * var = getVar(var_name, comp);
1562 0 : _c_fe_problem.needsPreviousNewtonIteration(true);
1563 0 : if (!var)
1564 : {
1565 0 : _default_gradient.resize(_coupleable_max_qps);
1566 0 : return _default_gradient;
1567 : }
1568 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1569 :
1570 0 : if (!_coupleable_neighbor)
1571 0 : return var->gradSlnPreviousNL();
1572 0 : return var->gradSlnPreviousNLNeighbor();
1573 : }
1574 :
1575 : const VariableGradient &
1576 25 : Coupleable::coupledGradientDot(const std::string & var_name, unsigned int comp) const
1577 : {
1578 25 : const auto * var = getVar(var_name, comp);
1579 25 : if (!var)
1580 : {
1581 0 : _default_gradient.resize(_coupleable_max_qps);
1582 0 : return _default_gradient;
1583 : }
1584 25 : checkFuncType(var_name, VarType::GradientDot, FuncAge::Curr);
1585 :
1586 25 : if (!_coupleable_neighbor)
1587 25 : return var->gradSlnDot();
1588 0 : return var->gradSlnNeighborDot();
1589 : }
1590 :
1591 : const VariableGradient &
1592 0 : Coupleable::coupledGradientDotDot(const std::string & var_name, unsigned int comp) const
1593 : {
1594 0 : const auto * var = getVar(var_name, comp);
1595 0 : if (!var)
1596 : {
1597 0 : _default_gradient.resize(_coupleable_max_qps);
1598 0 : return _default_gradient;
1599 : }
1600 0 : checkFuncType(var_name, VarType::GradientDot, FuncAge::Curr);
1601 :
1602 0 : if (!_coupleable_neighbor)
1603 0 : return var->gradSlnDotDot();
1604 0 : return var->gradSlnNeighborDotDot();
1605 : }
1606 :
1607 : const VectorVariableGradient &
1608 14 : Coupleable::coupledVectorGradient(const std::string & var_name, unsigned int comp) const
1609 : {
1610 14 : const auto * var = getVectorVar(var_name, comp);
1611 14 : if (!var)
1612 : {
1613 0 : _default_vector_gradient.resize(_coupleable_max_qps);
1614 0 : return _default_vector_gradient;
1615 : }
1616 14 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1617 :
1618 14 : if (!_coupleable_neighbor)
1619 14 : return (_c_is_implicit) ? var->gradSln() : var->gradSlnOld();
1620 0 : return (_c_is_implicit) ? var->gradSlnNeighbor() : var->gradSlnOldNeighbor();
1621 : }
1622 :
1623 : const VectorVariableGradient &
1624 14 : Coupleable::coupledVectorGradientOld(const std::string & var_name, unsigned int comp) const
1625 : {
1626 14 : const auto * var = getVectorVar(var_name, comp);
1627 14 : if (!var)
1628 : {
1629 0 : _default_vector_gradient.resize(_coupleable_max_qps);
1630 0 : return _default_vector_gradient;
1631 : }
1632 14 : checkFuncType(var_name, VarType::Gradient, FuncAge::Old);
1633 :
1634 14 : if (!_coupleable_neighbor)
1635 14 : return (_c_is_implicit) ? var->gradSlnOld() : var->gradSlnOlder();
1636 0 : return (_c_is_implicit) ? var->gradSlnOldNeighbor() : var->gradSlnOlderNeighbor();
1637 : }
1638 :
1639 : const VectorVariableGradient &
1640 14 : Coupleable::coupledVectorGradientOlder(const std::string & var_name, unsigned int comp) const
1641 : {
1642 14 : const auto * var = getVectorVar(var_name, comp);
1643 14 : if (!var)
1644 : {
1645 0 : _default_vector_gradient.resize(_coupleable_max_qps);
1646 0 : return _default_vector_gradient;
1647 : }
1648 14 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1649 :
1650 14 : if (!_coupleable_neighbor)
1651 14 : return var->gradSlnOlder();
1652 0 : return var->gradSlnOlderNeighbor();
1653 : }
1654 :
1655 : const ArrayVariableGradient &
1656 906 : Coupleable::coupledArrayGradient(const std::string & var_name, unsigned int comp) const
1657 : {
1658 906 : const auto * var = getArrayVar(var_name, comp);
1659 906 : if (!var)
1660 0 : return _default_array_gradient;
1661 906 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1662 :
1663 906 : if (!_coupleable_neighbor)
1664 906 : return (_c_is_implicit) ? var->gradSln() : var->gradSlnOld();
1665 0 : return (_c_is_implicit) ? var->gradSlnNeighbor() : var->gradSlnOldNeighbor();
1666 : }
1667 :
1668 : const ArrayVariableGradient &
1669 0 : Coupleable::coupledArrayGradientOld(const std::string & var_name, unsigned int comp) const
1670 : {
1671 0 : const auto * var = getArrayVar(var_name, comp);
1672 0 : if (!var)
1673 0 : return _default_array_gradient;
1674 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Old);
1675 :
1676 0 : if (!_coupleable_neighbor)
1677 0 : return (_c_is_implicit) ? var->gradSlnOld() : var->gradSlnOlder();
1678 0 : return (_c_is_implicit) ? var->gradSlnOldNeighbor() : var->gradSlnOlderNeighbor();
1679 : }
1680 :
1681 : const ArrayVariableGradient &
1682 0 : Coupleable::coupledArrayGradientOlder(const std::string & var_name, unsigned int comp) const
1683 : {
1684 0 : const auto * var = getArrayVar(var_name, comp);
1685 0 : if (!var)
1686 0 : return _default_array_gradient;
1687 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1688 :
1689 0 : if (!_coupleable_neighbor)
1690 0 : return var->gradSlnOlder();
1691 0 : return var->gradSlnOlderNeighbor();
1692 : }
1693 :
1694 : const ArrayVariableGradient &
1695 14 : Coupleable::coupledArrayGradientDot(const std::string & var_name, unsigned int comp) const
1696 : {
1697 14 : const auto * const var = getArrayVar(var_name, comp);
1698 14 : if (!var)
1699 0 : return _default_array_gradient;
1700 14 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1701 :
1702 14 : if (!_coupleable_neighbor)
1703 14 : return var->gradSlnDot();
1704 0 : return var->gradSlnNeighborDot();
1705 : }
1706 :
1707 : const VectorVariableCurl &
1708 0 : Coupleable::coupledCurl(const std::string & var_name, unsigned int comp) const
1709 : {
1710 0 : const auto * var = getVectorVar(var_name, comp);
1711 0 : if (!var)
1712 : {
1713 0 : _default_vector_curl.resize(_coupleable_max_qps);
1714 0 : return _default_vector_curl;
1715 : }
1716 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1717 :
1718 0 : if (!_coupleable_neighbor)
1719 0 : return (_c_is_implicit) ? var->curlSln() : var->curlSlnOld();
1720 0 : return (_c_is_implicit) ? var->curlSlnNeighbor() : var->curlSlnOldNeighbor();
1721 : }
1722 :
1723 : const VectorVariableCurl &
1724 0 : Coupleable::coupledCurlOld(const std::string & var_name, unsigned int comp) const
1725 : {
1726 0 : const auto * var = getVectorVar(var_name, comp);
1727 0 : if (!var)
1728 : {
1729 0 : _default_vector_curl.resize(_coupleable_max_qps);
1730 0 : return _default_vector_curl;
1731 : }
1732 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Old);
1733 :
1734 0 : if (!_coupleable_neighbor)
1735 0 : return (_c_is_implicit) ? var->curlSlnOld() : var->curlSlnOlder();
1736 0 : return (_c_is_implicit) ? var->curlSlnOldNeighbor() : var->curlSlnOlderNeighbor();
1737 : }
1738 :
1739 : const VectorVariableCurl &
1740 0 : Coupleable::coupledCurlOlder(const std::string & var_name, unsigned int comp) const
1741 : {
1742 0 : const auto * var = getVectorVar(var_name, comp);
1743 0 : if (!var)
1744 : {
1745 0 : _default_vector_curl.resize(_coupleable_max_qps);
1746 0 : return _default_vector_curl;
1747 : }
1748 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1749 :
1750 0 : if (!_coupleable_neighbor)
1751 0 : return var->curlSlnOlder();
1752 0 : return var->curlSlnOlderNeighbor();
1753 : }
1754 :
1755 : const ADVectorVariableCurl &
1756 14 : Coupleable::adCoupledCurl(const std::string & var_name, unsigned int comp) const
1757 : {
1758 14 : const auto * var = getVectorVar(var_name, comp);
1759 :
1760 14 : if (!var)
1761 0 : return getADDefaultCurl();
1762 14 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1763 :
1764 14 : if (!_c_is_implicit)
1765 0 : mooseError("Not implemented");
1766 :
1767 14 : if (!_coupleable_neighbor)
1768 14 : return var->adCurlSln();
1769 0 : return var->adCurlSlnNeighbor();
1770 : }
1771 :
1772 : const VectorVariableDivergence &
1773 168 : Coupleable::coupledDiv(const std::string & var_name, unsigned int comp) const
1774 : {
1775 168 : const auto * var = getVectorVar(var_name, comp);
1776 168 : if (!var)
1777 : {
1778 0 : _default_div.resize(_coupleable_max_qps);
1779 0 : return _default_div;
1780 : }
1781 168 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
1782 :
1783 168 : if (!_coupleable_neighbor)
1784 168 : return (_c_is_implicit) ? var->divSln() : var->divSlnOld();
1785 0 : return (_c_is_implicit) ? var->divSlnNeighbor() : var->divSlnOldNeighbor();
1786 : }
1787 :
1788 : const VectorVariableDivergence &
1789 0 : Coupleable::coupledDivOld(const std::string & var_name, unsigned int comp) const
1790 : {
1791 0 : const auto * var = getVectorVar(var_name, comp);
1792 0 : if (!var)
1793 : {
1794 0 : _default_div.resize(_coupleable_max_qps);
1795 0 : return _default_div;
1796 : }
1797 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Old);
1798 :
1799 0 : if (!_coupleable_neighbor)
1800 0 : return (_c_is_implicit) ? var->divSlnOld() : var->divSlnOlder();
1801 0 : return (_c_is_implicit) ? var->divSlnOldNeighbor() : var->divSlnOlderNeighbor();
1802 : }
1803 :
1804 : const VectorVariableDivergence &
1805 0 : Coupleable::coupledDivOlder(const std::string & var_name, unsigned int comp) const
1806 : {
1807 0 : const auto * var = getVectorVar(var_name, comp);
1808 0 : if (!var)
1809 : {
1810 0 : _default_div.resize(_coupleable_max_qps);
1811 0 : return _default_div;
1812 : }
1813 0 : checkFuncType(var_name, VarType::Gradient, FuncAge::Older);
1814 :
1815 0 : if (!_coupleable_neighbor)
1816 0 : return var->divSlnOlder();
1817 0 : return var->divSlnOlderNeighbor();
1818 : }
1819 :
1820 : const VariableSecond &
1821 52 : Coupleable::coupledSecond(const std::string & var_name, unsigned int comp) const
1822 : {
1823 52 : const auto * var = getVar(var_name, comp);
1824 52 : if (!var)
1825 : {
1826 38 : _default_second.resize(_coupleable_max_qps);
1827 38 : return _default_second;
1828 : }
1829 14 : checkFuncType(var_name, VarType::Second, FuncAge::Curr);
1830 :
1831 14 : if (!_coupleable_neighbor)
1832 14 : return (_c_is_implicit) ? var->secondSln() : var->secondSlnOlder();
1833 0 : return (_c_is_implicit) ? var->secondSlnNeighbor() : var->secondSlnOlderNeighbor();
1834 : }
1835 :
1836 : const VariableSecond &
1837 0 : Coupleable::coupledSecondOld(const std::string & var_name, unsigned int comp) const
1838 : {
1839 0 : const auto * var = getVar(var_name, comp);
1840 0 : if (!var)
1841 : {
1842 0 : _default_second.resize(_coupleable_max_qps);
1843 0 : return _default_second;
1844 : }
1845 0 : checkFuncType(var_name, VarType::Second, FuncAge::Old);
1846 :
1847 0 : if (!_coupleable_neighbor)
1848 0 : return (_c_is_implicit) ? var->secondSlnOld() : var->secondSlnOlder();
1849 0 : return (_c_is_implicit) ? var->secondSlnOldNeighbor() : var->secondSlnOlderNeighbor();
1850 : }
1851 :
1852 : const VariableSecond &
1853 0 : Coupleable::coupledSecondOlder(const std::string & var_name, unsigned int comp) const
1854 : {
1855 0 : const auto * var = getVar(var_name, comp);
1856 0 : if (!var)
1857 : {
1858 0 : _default_second.resize(_coupleable_max_qps);
1859 0 : return _default_second;
1860 : }
1861 0 : checkFuncType(var_name, VarType::Second, FuncAge::Older);
1862 :
1863 0 : if (!_coupleable_neighbor)
1864 0 : return var->secondSlnOlder();
1865 0 : return var->secondSlnOlderNeighbor();
1866 : }
1867 :
1868 : const VariableSecond &
1869 0 : Coupleable::coupledSecondPreviousNL(const std::string & var_name, unsigned int comp) const
1870 : {
1871 0 : const auto * var = getVar(var_name, comp);
1872 0 : _c_fe_problem.needsPreviousNewtonIteration(true);
1873 0 : if (!var)
1874 : {
1875 0 : _default_second.resize(_coupleable_max_qps);
1876 0 : return _default_second;
1877 : }
1878 0 : checkFuncType(var_name, VarType::Second, FuncAge::Curr);
1879 :
1880 0 : if (!_coupleable_neighbor)
1881 0 : return var->secondSlnPreviousNL();
1882 0 : return var->secondSlnPreviousNLNeighbor();
1883 : }
1884 :
1885 : template <typename T>
1886 : const T &
1887 46 : Coupleable::coupledNodalValue(const std::string & var_name, unsigned int comp) const
1888 : {
1889 46 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
1890 46 : if (!var)
1891 0 : return getDefaultNodalValue<T>(var_name, comp);
1892 46 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
1893 :
1894 46 : if (!var->isNodal())
1895 4 : mooseError(_c_name,
1896 : ": Trying to get nodal values of variable '",
1897 4 : var->name(),
1898 : "', but it is not nodal.");
1899 :
1900 42 : if (!_coupleable_neighbor)
1901 42 : return (_c_is_implicit) ? var->nodalValue() : var->nodalValueOld();
1902 0 : return (_c_is_implicit) ? var->nodalValueNeighbor() : var->nodalValueOldNeighbor();
1903 : }
1904 :
1905 : template <typename T>
1906 : const T &
1907 4 : Coupleable::coupledNodalValueOld(const std::string & var_name, unsigned int comp) const
1908 : {
1909 4 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
1910 4 : if (!var)
1911 0 : return getDefaultNodalValue<T>(var_name, comp);
1912 4 : checkFuncType(var_name, VarType::Ignore, FuncAge::Old);
1913 :
1914 4 : if (!var->isNodal())
1915 4 : mooseError(_c_name,
1916 : ": Trying to get old nodal values of variable '",
1917 4 : var->name(),
1918 : "', but it is not nodal.");
1919 :
1920 0 : if (!_coupleable_neighbor)
1921 0 : return (_c_is_implicit) ? var->nodalValueOld() : var->nodalValueOlder();
1922 0 : return (_c_is_implicit) ? var->nodalValueOldNeighbor() : var->nodalValueOlderNeighbor();
1923 : }
1924 :
1925 : template <typename T>
1926 : const T &
1927 4 : Coupleable::coupledNodalValueOlder(const std::string & var_name, unsigned int comp) const
1928 : {
1929 4 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
1930 4 : if (!var)
1931 0 : return getDefaultNodalValue<T>(var_name, comp);
1932 4 : checkFuncType(var_name, VarType::Ignore, FuncAge::Older);
1933 :
1934 4 : if (!var->isNodal())
1935 4 : mooseError(_c_name,
1936 : ": Trying to get older nodal values of variable '",
1937 4 : var->name(),
1938 : "', but it is not nodal.");
1939 :
1940 0 : if (!_coupleable_neighbor)
1941 0 : return var->nodalValueOlder();
1942 0 : return var->nodalValueOlderNeighbor();
1943 : }
1944 :
1945 : template <typename T>
1946 : const T &
1947 0 : Coupleable::coupledNodalValuePreviousNL(const std::string & var_name, unsigned int comp) const
1948 : {
1949 0 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
1950 0 : if (!var)
1951 0 : return getDefaultNodalValue<T>(var_name, comp);
1952 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
1953 :
1954 0 : _c_fe_problem.needsPreviousNewtonIteration(true);
1955 :
1956 0 : if (!_coupleable_neighbor)
1957 0 : return var->nodalValuePreviousNL();
1958 0 : return var->nodalValuePreviousNLNeighbor();
1959 : }
1960 :
1961 : template <typename T>
1962 : const T &
1963 0 : Coupleable::coupledNodalDot(const std::string & var_name, unsigned int comp) const
1964 : {
1965 0 : static const T zero = 0;
1966 0 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
1967 0 : if (!var)
1968 0 : return zero;
1969 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1970 :
1971 0 : if (!_coupleable_neighbor)
1972 0 : return var->nodalValueDot();
1973 0 : mooseError("Neighbor version not implemented");
1974 : }
1975 :
1976 : const VariableValue &
1977 0 : Coupleable::coupledNodalDotDot(const std::string & var_name, unsigned int comp) const
1978 : {
1979 0 : const auto * var = getVar(var_name, comp);
1980 0 : if (!var)
1981 : {
1982 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1983 0 : return _default_value_zero;
1984 : }
1985 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
1986 :
1987 0 : if (!_coupleable_neighbor)
1988 0 : return var->dofValuesDotDot();
1989 0 : return var->dofValuesDotDotNeighbor();
1990 : }
1991 :
1992 : const VariableValue &
1993 0 : Coupleable::coupledNodalDotOld(const std::string & var_name, unsigned int comp) const
1994 : {
1995 0 : const auto * var = getVar(var_name, comp);
1996 0 : if (!var)
1997 : {
1998 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
1999 0 : return _default_value_zero;
2000 : }
2001 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
2002 :
2003 0 : if (!_coupleable_neighbor)
2004 0 : return var->dofValuesDotOld();
2005 0 : return var->dofValuesDotOldNeighbor();
2006 : }
2007 :
2008 : const VariableValue &
2009 0 : Coupleable::coupledNodalDotDotOld(const std::string & var_name, unsigned int comp) const
2010 : {
2011 0 : const auto * var = getVar(var_name, comp);
2012 0 : if (!var)
2013 : {
2014 0 : _default_value_zero.resize(_coupleable_max_qps, 0);
2015 0 : return _default_value_zero;
2016 : }
2017 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Old);
2018 :
2019 0 : if (!_coupleable_neighbor)
2020 0 : return var->dofValuesDotDotOld();
2021 0 : return var->dofValuesDotDotOldNeighbor();
2022 : }
2023 :
2024 : const VariableValue &
2025 170 : Coupleable::coupledDofValues(const std::string & var_name, unsigned int comp) const
2026 : {
2027 170 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2028 170 : if (!var)
2029 0 : return *getDefaultValue(var_name, comp);
2030 170 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2031 :
2032 170 : if (!_coupleable_neighbor)
2033 156 : return (_c_is_implicit) ? var->dofValues() : var->dofValuesOld();
2034 14 : return (_c_is_implicit) ? var->dofValuesNeighbor() : var->dofValuesOldNeighbor();
2035 : }
2036 :
2037 : std::vector<const VariableValue *>
2038 36 : Coupleable::coupledAllDofValues(const std::string & var_name) const
2039 : {
2040 72 : auto func = [this, &var_name](unsigned int comp) { return &coupledDofValues(var_name, comp); };
2041 72 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2042 : }
2043 :
2044 : const VariableValue &
2045 14 : Coupleable::coupledDofValuesOld(const std::string & var_name, unsigned int comp) const
2046 : {
2047 14 : const auto * var = getVar(var_name, comp);
2048 14 : if (!var)
2049 0 : return *getDefaultValue(var_name, comp);
2050 14 : checkFuncType(var_name, VarType::Ignore, FuncAge::Old);
2051 :
2052 14 : if (!_coupleable_neighbor)
2053 14 : return (_c_is_implicit) ? var->dofValuesOld() : var->dofValuesOlder();
2054 0 : return (_c_is_implicit) ? var->dofValuesOldNeighbor() : var->dofValuesOlderNeighbor();
2055 : }
2056 :
2057 : std::vector<const VariableValue *>
2058 0 : Coupleable::coupledAllDofValuesOld(const std::string & var_name) const
2059 : {
2060 0 : auto func = [this, &var_name](unsigned int comp) { return &coupledDofValuesOld(var_name, comp); };
2061 0 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2062 : }
2063 :
2064 : const VariableValue &
2065 0 : Coupleable::coupledDofValuesOlder(const std::string & var_name, unsigned int comp) const
2066 : {
2067 0 : const auto * var = getVar(var_name, comp);
2068 0 : if (!var)
2069 0 : return *getDefaultValue(var_name, comp);
2070 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Older);
2071 :
2072 0 : if (!_coupleable_neighbor)
2073 0 : return var->dofValuesOlder();
2074 0 : return var->dofValuesOlderNeighbor();
2075 : }
2076 :
2077 : std::vector<const VariableValue *>
2078 0 : Coupleable::coupledAllDofValuesOlder(const std::string & var_name) const
2079 : {
2080 0 : auto func = [this, &var_name](unsigned int comp)
2081 0 : { return &coupledDofValuesOlder(var_name, comp); };
2082 0 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2083 : }
2084 :
2085 : const ArrayVariableValue &
2086 0 : Coupleable::coupledArrayDofValues(const std::string & var_name, unsigned int comp) const
2087 : {
2088 0 : const auto * var = getArrayVar(var_name, comp);
2089 0 : if (!var)
2090 0 : return *getDefaultArrayValue(var_name);
2091 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2092 :
2093 0 : if (!_coupleable_neighbor)
2094 0 : return (_c_is_implicit) ? var->dofValues() : var->dofValuesOld();
2095 0 : return (_c_is_implicit) ? var->dofValuesNeighbor() : var->dofValuesOldNeighbor();
2096 : }
2097 :
2098 : const ADVariableValue &
2099 0 : Coupleable::adCoupledDofValues(const std::string & var_name, unsigned int comp) const
2100 : {
2101 0 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2102 :
2103 0 : if (!var)
2104 0 : return *getADDefaultValue(var_name);
2105 0 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2106 :
2107 0 : if (!_c_is_implicit)
2108 0 : mooseError("Not implemented");
2109 :
2110 0 : if (!_coupleable_neighbor)
2111 0 : return var->adDofValues();
2112 0 : return var->adDofValuesNeighbor();
2113 : }
2114 :
2115 : void
2116 25891 : Coupleable::validateExecutionerType(const std::string & name, const std::string & fn_name) const
2117 : {
2118 25891 : if (!_c_fe_problem.isTransient())
2119 8 : mooseError(_c_name,
2120 : ": Calling \"",
2121 : fn_name,
2122 : "\" on variable \"",
2123 : name,
2124 : "\" when using a \"Steady\" executioner is not allowed. This value is available "
2125 : "only in transient simulations.");
2126 25883 : }
2127 :
2128 : template <typename T>
2129 : const typename Moose::ADType<T>::type &
2130 23 : Coupleable::adCoupledNodalValue(const std::string & var_name, unsigned int comp) const
2131 : {
2132 23 : static const typename Moose::ADType<T>::type zero = 0;
2133 23 : if (!isCoupled(var_name))
2134 0 : return zero;
2135 :
2136 23 : if (!_c_nodal)
2137 0 : mooseError("The adCoupledNodalValue method should only be called for nodal computing objects");
2138 23 : if (_coupleable_neighbor)
2139 0 : mooseError(
2140 : "The adCoupledNodalValue method shouldn't be called for neighbor computing objects. I "
2141 : "don't even know what that would mean, although maybe someone could explain it to me.");
2142 23 : if (!_c_is_implicit)
2143 0 : mooseError("If you're going to use an explicit scheme, then use coupledNodalValue instead of "
2144 : "adCoupledNodalValue");
2145 :
2146 23 : const auto * var = getVarHelper<MooseVariableFE<T>>(var_name, comp);
2147 :
2148 23 : return var->adNodalValue();
2149 : }
2150 :
2151 : const ADVariableValue &
2152 3338 : Coupleable::adCoupledValue(const std::string & var_name, unsigned int comp) const
2153 : {
2154 3338 : const auto * const var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2155 :
2156 3338 : if (!var)
2157 147 : return *getADDefaultValue(var_name);
2158 3191 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2159 :
2160 3191 : if (!_c_is_implicit)
2161 0 : mooseError("Not implemented");
2162 :
2163 3191 : if (_c_nodal)
2164 131 : return var->adDofValues();
2165 :
2166 3060 : if (!_coupleable_neighbor)
2167 2233 : return var->adSln();
2168 827 : return var->adSlnNeighbor();
2169 : }
2170 :
2171 : const ADVariableValue &
2172 28 : Coupleable::adCoupledLowerValue(const std::string & var_name, unsigned int comp) const
2173 : {
2174 28 : auto var = getVarHelper<MooseVariableFE<Real>>(var_name, comp);
2175 :
2176 28 : if (!var)
2177 0 : return *getADDefaultValue(var_name);
2178 28 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2179 :
2180 28 : if (!_c_is_implicit)
2181 0 : mooseError("adCoupledLowerValue cannot be called in a coupleable neighbor object");
2182 :
2183 28 : if (_c_nodal)
2184 0 : return var->adDofValues();
2185 : else
2186 28 : return var->adSlnLower();
2187 : }
2188 :
2189 : const ADVariableGradient &
2190 198 : Coupleable::adCoupledGradient(const std::string & var_name, unsigned int comp) const
2191 : {
2192 198 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2193 :
2194 198 : if (!var)
2195 0 : return getADDefaultGradient();
2196 198 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
2197 :
2198 198 : if (!_c_is_implicit)
2199 0 : mooseError("Not implemented");
2200 :
2201 198 : if (!_coupleable_neighbor)
2202 156 : return var->adGradSln();
2203 42 : return var->adGradSlnNeighbor();
2204 : }
2205 :
2206 : const ADVariableGradient &
2207 0 : Coupleable::adCoupledGradientDot(const std::string & var_name, unsigned int comp) const
2208 : {
2209 0 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2210 :
2211 0 : if (!var)
2212 0 : return getADDefaultGradient();
2213 0 : checkFuncType(var_name, VarType::GradientDot, FuncAge::Curr);
2214 :
2215 0 : if (!_c_is_implicit)
2216 0 : mooseError("Not implemented");
2217 :
2218 0 : if (!_coupleable_neighbor)
2219 0 : return var->adGradSlnDot();
2220 0 : return var->adGradSlnNeighborDot();
2221 : }
2222 :
2223 : const ADVariableSecond &
2224 0 : Coupleable::adCoupledSecond(const std::string & var_name, unsigned int comp) const
2225 : {
2226 0 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2227 :
2228 0 : if (!var)
2229 0 : return getADDefaultSecond();
2230 0 : checkFuncType(var_name, VarType::Second, FuncAge::Curr);
2231 :
2232 0 : if (!_c_is_implicit)
2233 0 : mooseError("Not implemented");
2234 :
2235 0 : if (!_coupleable_neighbor)
2236 0 : return var->adSecondSln();
2237 : else
2238 0 : return var->adSecondSlnNeighbor();
2239 : }
2240 :
2241 : const ADVectorVariableSecond &
2242 0 : adCoupledVectorSecond(const std::string & /*var_name*/, unsigned int /*comp = 0*/)
2243 : {
2244 0 : mooseError("Automatic differentiation using second derivatives of vector variables is not "
2245 : "implemented.");
2246 : }
2247 :
2248 : const ADVariableValue &
2249 84 : Coupleable::adCoupledDot(const std::string & var_name, unsigned int comp) const
2250 : {
2251 84 : const auto * var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2252 :
2253 84 : if (!var)
2254 0 : return *getADDefaultValue(var_name);
2255 84 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
2256 :
2257 84 : if (!_coupleable_neighbor)
2258 : {
2259 70 : if (_c_nodal)
2260 14 : return var->adDofValuesDot();
2261 56 : return var->adUDot();
2262 : }
2263 : else
2264 : {
2265 14 : if (_c_nodal)
2266 0 : mooseError("AD neighbor nodal dof dot not implemented");
2267 14 : return var->adUDotNeighbor();
2268 : }
2269 : }
2270 :
2271 : const ADVariableValue &
2272 70 : Coupleable::adCoupledDotDot(const std::string & var_name, unsigned int comp) const
2273 : {
2274 70 : const auto * const var = getVarHelper<MooseVariableField<Real>>(var_name, comp);
2275 :
2276 70 : if (!var)
2277 0 : return *getADDefaultValue(var_name);
2278 70 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
2279 :
2280 70 : if (_c_nodal)
2281 0 : mooseError("Not implemented");
2282 :
2283 70 : if (!_coupleable_neighbor)
2284 56 : return var->adUDotDot();
2285 14 : return var->adUDotDotNeighbor();
2286 : }
2287 :
2288 : const ADVectorVariableValue &
2289 0 : Coupleable::adCoupledVectorDot(const std::string & var_name, unsigned int comp) const
2290 : {
2291 0 : const auto * var = getVectorVar(var_name, comp);
2292 0 : if (!var)
2293 0 : return *getADDefaultVectorValue(var_name);
2294 0 : checkFuncType(var_name, VarType::Dot, FuncAge::Curr);
2295 :
2296 0 : if (_c_nodal)
2297 0 : mooseError("Not implemented");
2298 :
2299 0 : if (!_coupleable_neighbor)
2300 0 : return var->adUDot();
2301 0 : return var->adUDotNeighbor();
2302 : }
2303 :
2304 : const ADVectorVariableValue &
2305 165 : Coupleable::adCoupledVectorValue(const std::string & var_name, unsigned int comp) const
2306 : {
2307 165 : const auto * var = getVectorVar(var_name, comp);
2308 165 : if (!var)
2309 28 : return *getADDefaultVectorValue(var_name);
2310 137 : checkFuncType(var_name, VarType::Ignore, FuncAge::Curr);
2311 :
2312 137 : if (_c_nodal)
2313 0 : mooseError("Not implemented");
2314 137 : if (!_c_is_implicit)
2315 0 : mooseError("Not implemented");
2316 :
2317 137 : if (!_coupleable_neighbor)
2318 109 : return var->adSln();
2319 28 : return var->adSlnNeighbor();
2320 : }
2321 :
2322 : const ADVectorVariableGradient &
2323 42 : Coupleable::adCoupledVectorGradient(const std::string & var_name, unsigned int comp) const
2324 : {
2325 42 : const auto * var = getVectorVar(var_name, comp);
2326 42 : if (!var)
2327 14 : return getADDefaultVectorGradient();
2328 28 : checkFuncType(var_name, VarType::Gradient, FuncAge::Curr);
2329 :
2330 28 : if (!_c_is_implicit)
2331 0 : mooseError("Not implemented");
2332 :
2333 28 : if (!_coupleable_neighbor)
2334 28 : return var->adGradSln();
2335 0 : return var->adGradSlnNeighbor();
2336 : }
2337 :
2338 : const ADVariableValue *
2339 147 : Coupleable::getADDefaultValue(const std::string & var_name) const
2340 : {
2341 147 : auto default_value_it = _ad_default_value.find(var_name);
2342 147 : if (default_value_it == _ad_default_value.end())
2343 : {
2344 147 : auto value = std::make_unique<ADVariableValue>(_coupleable_max_qps,
2345 147 : _c_parameters.defaultCoupledValue(var_name));
2346 147 : default_value_it = _ad_default_value.insert(std::make_pair(var_name, std::move(value))).first;
2347 147 : }
2348 :
2349 294 : return default_value_it->second.get();
2350 : }
2351 :
2352 : const ADVectorVariableValue *
2353 28 : Coupleable::getADDefaultVectorValue(const std::string & var_name) const
2354 : {
2355 28 : auto default_value_it = _ad_default_vector_value.find(var_name);
2356 28 : if (default_value_it == _ad_default_vector_value.end())
2357 : {
2358 28 : RealVectorValue default_vec;
2359 98 : for (unsigned int i = 0; i < _c_parameters.numberDefaultCoupledValues(var_name); ++i)
2360 70 : default_vec(i) = _c_parameters.defaultCoupledValue(var_name, i);
2361 28 : auto value = std::make_unique<ADVectorVariableValue>(_coupleable_max_qps, default_vec);
2362 28 : default_value_it =
2363 28 : _ad_default_vector_value.insert(std::make_pair(var_name, std::move(value))).first;
2364 28 : }
2365 :
2366 56 : return default_value_it->second.get();
2367 : }
2368 :
2369 : const ADVariableGradient &
2370 0 : Coupleable::getADDefaultGradient() const
2371 : {
2372 0 : _ad_default_gradient.resize(_coupleable_max_qps);
2373 0 : return _ad_default_gradient;
2374 : }
2375 :
2376 : const ADVectorVariableGradient &
2377 14 : Coupleable::getADDefaultVectorGradient() const
2378 : {
2379 14 : _ad_default_vector_gradient.resize(_coupleable_max_qps);
2380 14 : return _ad_default_vector_gradient;
2381 : }
2382 :
2383 : const ADVariableSecond &
2384 0 : Coupleable::getADDefaultSecond() const
2385 : {
2386 0 : _ad_default_second.resize(_coupleable_max_qps);
2387 0 : return _ad_default_second;
2388 : }
2389 :
2390 : const ADVectorVariableCurl &
2391 0 : Coupleable::getADDefaultCurl() const
2392 : {
2393 0 : _ad_default_curl.resize(_coupleable_max_qps);
2394 0 : return _ad_default_curl;
2395 : }
2396 :
2397 : const ADVariableValue &
2398 0 : Coupleable::adZeroValue() const
2399 : {
2400 0 : mooseDeprecated("Method adZeroValue() is deprecated. Use '_ad_zero' instead.");
2401 0 : return _ad_zero;
2402 : }
2403 :
2404 : const ADVariableGradient &
2405 0 : Coupleable::adZeroGradient() const
2406 : {
2407 0 : mooseDeprecated("Method adZeroGradient() is deprecated. Use '_ad_grad_zero' instead.");
2408 0 : return _ad_grad_zero;
2409 : }
2410 :
2411 : const ADVariableSecond &
2412 0 : Coupleable::adZeroSecond() const
2413 : {
2414 0 : mooseDeprecated("Method adZeroSecond() is deprecated. Use '_ad_second_zero' instead.");
2415 0 : return _ad_second_zero;
2416 : }
2417 :
2418 : template <>
2419 : const GenericVariableValue<false> &
2420 42 : Coupleable::genericZeroValue<false>()
2421 : {
2422 42 : return _zero;
2423 : }
2424 :
2425 : template <>
2426 : const GenericVariableValue<true> &
2427 84 : Coupleable::genericZeroValue<true>()
2428 : {
2429 84 : return _ad_zero;
2430 : }
2431 :
2432 : template <>
2433 : const GenericVariableGradient<false> &
2434 0 : Coupleable::genericZeroGradient<false>()
2435 : {
2436 0 : return _grad_zero;
2437 : }
2438 :
2439 : template <>
2440 : const GenericVariableGradient<true> &
2441 0 : Coupleable::genericZeroGradient<true>()
2442 : {
2443 0 : return _ad_grad_zero;
2444 : }
2445 :
2446 : template <>
2447 : const GenericVariableSecond<false> &
2448 0 : Coupleable::genericZeroSecond<false>()
2449 : {
2450 0 : return _second_zero;
2451 : }
2452 :
2453 : template <>
2454 : const GenericVariableSecond<true> &
2455 0 : Coupleable::genericZeroSecond<true>()
2456 : {
2457 0 : return _ad_second_zero;
2458 : }
2459 :
2460 : template <>
2461 : const GenericVariableGradient<false> &
2462 42 : Coupleable::coupledGenericGradient<false>(const std::string & var_name, unsigned int comp) const
2463 : {
2464 42 : return coupledGradient(var_name, comp);
2465 : }
2466 :
2467 : template <>
2468 : const GenericVariableGradient<true> &
2469 126 : Coupleable::coupledGenericGradient<true>(const std::string & var_name, unsigned int comp) const
2470 : {
2471 126 : return adCoupledGradient(var_name, comp);
2472 : }
2473 :
2474 : std::vector<unsigned int>
2475 67 : Coupleable::coupledIndices(const std::string & var_name) const
2476 : {
2477 109 : auto func = [this, &var_name](unsigned int comp) { return coupled(var_name, comp); };
2478 134 : return coupledVectorHelper<unsigned int>(var_name, func);
2479 : }
2480 :
2481 : VariableName
2482 5520 : Coupleable::coupledName(const std::string & var_name, unsigned int comp) const
2483 : {
2484 5520 : if (getFieldVar(var_name, comp))
2485 5516 : return getFieldVar(var_name, comp)->name();
2486 : // Detect if we are in the case where a constant was passed in lieu of a variable
2487 4 : else if (isCoupledConstant(var_name))
2488 4 : mooseError(_c_name,
2489 : ": a variable name was queried but a constant was passed for parameter '",
2490 : var_name,
2491 : "Either pass a true variable or contact a developer to shield the call to "
2492 : "'coupledName' with 'isCoupledConstant'");
2493 : else
2494 0 : mooseError(
2495 0 : _c_name, ": Variable '", var_name, "' does not exist, yet its coupled name is requested");
2496 : }
2497 :
2498 : std::vector<VariableName>
2499 0 : Coupleable::coupledNames(const std::string & var_name) const
2500 : {
2501 0 : auto func = [this, &var_name](unsigned int comp) { return coupledName(var_name, comp); };
2502 0 : return coupledVectorHelper<VariableName>(var_name, func);
2503 : }
2504 :
2505 : std::vector<const VariableValue *>
2506 1716 : Coupleable::coupledValues(const std::string & var_name) const
2507 : {
2508 1274 : auto func = [this, &var_name](unsigned int comp) { return &coupledValue(var_name, comp); };
2509 3432 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2510 : }
2511 :
2512 : std::vector<const VectorVariableValue *>
2513 42 : Coupleable::coupledVectorValues(const std::string & var_name) const
2514 : {
2515 56 : auto func = [this, &var_name](unsigned int comp) { return &coupledVectorValue(var_name, comp); };
2516 84 : return coupledVectorHelper<const VectorVariableValue *>(var_name, func);
2517 : }
2518 :
2519 : template <>
2520 : std::vector<const GenericVariableValue<false> *>
2521 92 : Coupleable::coupledGenericValues<false>(const std::string & var_name) const
2522 : {
2523 92 : return coupledValues(var_name);
2524 : }
2525 :
2526 : template <>
2527 : std::vector<const GenericVariableValue<true> *>
2528 66 : Coupleable::coupledGenericValues<true>(const std::string & var_name) const
2529 : {
2530 66 : return adCoupledValues(var_name);
2531 : }
2532 :
2533 : std::vector<const ADVariableValue *>
2534 66 : Coupleable::adCoupledValues(const std::string & var_name) const
2535 : {
2536 66 : auto func = [this, &var_name](unsigned int comp) { return &adCoupledValue(var_name, comp); };
2537 132 : return coupledVectorHelper<const ADVariableValue *>(var_name, func);
2538 : }
2539 :
2540 : std::vector<const ADVectorVariableValue *>
2541 0 : Coupleable::adCoupledVectorValues(const std::string & var_name) const
2542 : {
2543 0 : auto func = [this, &var_name](unsigned int comp)
2544 0 : { return &adCoupledVectorValue(var_name, comp); };
2545 0 : return coupledVectorHelper<const ADVectorVariableValue *>(var_name, func);
2546 : }
2547 :
2548 : std::vector<const VariableValue *>
2549 14 : Coupleable::coupledVectorTagValues(const std::string & var_names, TagID tag) const
2550 : {
2551 14 : auto func = [this, &var_names, &tag](unsigned int comp)
2552 14 : { return &coupledVectorTagValue(var_names, tag, comp); };
2553 28 : return coupledVectorHelper<const VariableValue *>(var_names, func);
2554 : }
2555 :
2556 : std::vector<const VariableValue *>
2557 14 : Coupleable::coupledVectorTagValues(const std::string & var_names,
2558 : const std::string & tag_name) const
2559 : {
2560 14 : if (!_c_parameters.isParamValid(tag_name))
2561 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2562 :
2563 14 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2564 14 : if (!_c_fe_problem.vectorTagExists(tagname))
2565 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
2566 :
2567 14 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
2568 28 : return coupledVectorTagValues(var_names, tag);
2569 14 : }
2570 :
2571 : std::vector<const ArrayVariableValue *>
2572 0 : Coupleable::coupledVectorTagArrayValues(const std::string & var_names, TagID tag) const
2573 : {
2574 0 : auto func = [this, &var_names, &tag](unsigned int index)
2575 0 : { return &coupledVectorTagArrayValue(var_names, tag, index); };
2576 0 : return coupledVectorHelper<const ArrayVariableValue *>(var_names, func);
2577 : }
2578 :
2579 : std::vector<const ArrayVariableValue *>
2580 0 : Coupleable::coupledVectorTagArrayValues(const std::string & var_names,
2581 : const std::string & tag_name) const
2582 : {
2583 0 : if (!_c_parameters.isParamValid(tag_name))
2584 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2585 :
2586 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2587 0 : if (!_c_fe_problem.vectorTagExists(tagname))
2588 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
2589 :
2590 0 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
2591 0 : return coupledVectorTagArrayValues(var_names, tag);
2592 0 : }
2593 :
2594 : std::vector<const VariableGradient *>
2595 0 : Coupleable::coupledVectorTagGradients(const std::string & var_names, TagID tag) const
2596 : {
2597 0 : auto func = [this, &var_names, &tag](unsigned int index)
2598 0 : { return &coupledVectorTagGradient(var_names, tag, index); };
2599 0 : return coupledVectorHelper<const VariableGradient *>(var_names, func);
2600 : }
2601 :
2602 : std::vector<const VariableGradient *>
2603 0 : Coupleable::coupledVectorTagGradients(const std::string & var_names,
2604 : const std::string & tag_name) const
2605 : {
2606 0 : if (!_c_parameters.isParamValid(tag_name))
2607 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2608 :
2609 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2610 0 : if (!_c_fe_problem.vectorTagExists(tagname))
2611 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
2612 :
2613 0 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
2614 0 : return coupledVectorTagGradients(var_names, tag);
2615 0 : }
2616 :
2617 : std::vector<const ArrayVariableGradient *>
2618 0 : Coupleable::coupledVectorTagArrayGradients(const std::string & var_names, TagID tag) const
2619 : {
2620 0 : auto func = [this, &var_names, &tag](unsigned int index)
2621 0 : { return &coupledVectorTagArrayGradient(var_names, tag, index); };
2622 0 : return coupledVectorHelper<const ArrayVariableGradient *>(var_names, func);
2623 : }
2624 :
2625 : std::vector<const ArrayVariableGradient *>
2626 0 : Coupleable::coupledVectorTagArrayGradients(const std::string & var_names,
2627 : const std::string & tag_name) const
2628 : {
2629 0 : if (!_c_parameters.isParamValid(tag_name))
2630 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2631 :
2632 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2633 0 : if (!_c_fe_problem.vectorTagExists(tagname))
2634 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
2635 :
2636 0 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
2637 0 : return coupledVectorTagArrayGradients(var_names, tag);
2638 0 : }
2639 :
2640 : std::vector<const VariableValue *>
2641 0 : Coupleable::coupledVectorTagDofValues(const std::string & var_names, TagID tag) const
2642 : {
2643 0 : auto func = [this, &var_names, &tag](unsigned int comp)
2644 0 : { return &coupledVectorTagDofValue(var_names, tag, comp); };
2645 0 : return coupledVectorHelper<const VariableValue *>(var_names, func);
2646 : }
2647 :
2648 : std::vector<const VariableValue *>
2649 0 : Coupleable::coupledVectorTagDofValues(const std::string & var_names,
2650 : const std::string & tag_name) const
2651 : {
2652 0 : if (!_c_parameters.isParamValid(tag_name))
2653 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2654 :
2655 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2656 0 : if (!_c_fe_problem.vectorTagExists(tagname))
2657 0 : mooseError("Tagged vector with tag name '", tagname, "' does not exist");
2658 :
2659 0 : TagID tag = _c_fe_problem.getVectorTagID(tagname);
2660 0 : return coupledVectorTagDofValues(var_names, tag);
2661 0 : }
2662 :
2663 : std::vector<const VariableValue *>
2664 0 : Coupleable::coupledMatrixTagValues(const std::string & var_names, TagID tag) const
2665 : {
2666 0 : auto func = [this, &var_names, &tag](unsigned int comp)
2667 0 : { return &coupledMatrixTagValue(var_names, tag, comp); };
2668 0 : return coupledVectorHelper<const VariableValue *>(var_names, func);
2669 : }
2670 :
2671 : std::vector<const VariableValue *>
2672 0 : Coupleable::coupledMatrixTagValues(const std::string & var_names,
2673 : const std::string & tag_name) const
2674 : {
2675 0 : if (!_c_parameters.isParamValid(tag_name))
2676 0 : mooseError("Tag name parameter '", tag_name, "' is invalid");
2677 :
2678 0 : TagName tagname = _c_parameters.get<TagName>(tag_name);
2679 0 : if (!_c_fe_problem.matrixTagExists(tagname))
2680 0 : mooseError("Matrix tag name '", tagname, "' does not exist");
2681 :
2682 0 : TagID tag = _c_fe_problem.getMatrixTagID(tagname);
2683 0 : return coupledMatrixTagValues(var_names, tag);
2684 0 : }
2685 :
2686 : std::vector<const VariableValue *>
2687 504 : Coupleable::coupledValuesOld(const std::string & var_name) const
2688 : {
2689 11844 : auto func = [this, &var_name](unsigned int comp) { return &coupledValueOld(var_name, comp); };
2690 1008 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2691 : }
2692 :
2693 : std::vector<const VariableValue *>
2694 504 : Coupleable::coupledValuesOlder(const std::string & var_name) const
2695 : {
2696 11844 : auto func = [this, &var_name](unsigned int comp) { return &coupledValueOlder(var_name, comp); };
2697 1008 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2698 : }
2699 :
2700 : std::vector<const VectorVariableValue *>
2701 28 : Coupleable::coupledVectorValuesOld(const std::string & var_name) const
2702 : {
2703 56 : auto func = [this, &var_name](unsigned int comp)
2704 56 : { return &coupledVectorValueOld(var_name, comp); };
2705 56 : return coupledVectorHelper<const VectorVariableValue *>(var_name, func);
2706 : }
2707 :
2708 : std::vector<const VariableGradient *>
2709 0 : Coupleable::coupledGradients(const std::string & var_name) const
2710 : {
2711 0 : auto func = [this, &var_name](unsigned int comp) { return &coupledGradient(var_name, comp); };
2712 0 : return coupledVectorHelper<const VariableGradient *>(var_name, func);
2713 : }
2714 :
2715 : template <>
2716 : std::vector<const GenericVariableGradient<false> *>
2717 0 : Coupleable::coupledGenericGradients<false>(const std::string & var_name) const
2718 : {
2719 0 : return coupledGradients(var_name);
2720 : }
2721 :
2722 : template <>
2723 : std::vector<const GenericVariableGradient<true> *>
2724 0 : Coupleable::coupledGenericGradients<true>(const std::string & var_name) const
2725 : {
2726 0 : auto func = [this, &var_name](unsigned int comp) { return &adCoupledGradient(var_name, comp); };
2727 0 : return coupledVectorHelper<const GenericVariableGradient<true> *>(var_name, func);
2728 : }
2729 :
2730 : std::vector<const ADVariableGradient *>
2731 0 : Coupleable::adCoupledGradients(const std::string & var_name) const
2732 : {
2733 0 : auto func = [this, &var_name](unsigned int comp) { return &adCoupledGradient(var_name, comp); };
2734 0 : return coupledVectorHelper<const ADVariableGradient *>(var_name, func);
2735 : }
2736 :
2737 : std::vector<const VariableGradient *>
2738 0 : Coupleable::coupledGradientsOld(const std::string & var_name) const
2739 : {
2740 0 : auto func = [this, &var_name](unsigned int comp) { return &coupledGradientOld(var_name, comp); };
2741 0 : return coupledVectorHelper<const VariableGradient *>(var_name, func);
2742 : }
2743 :
2744 : std::vector<const VariableValue *>
2745 0 : Coupleable::coupledDots(const std::string & var_name) const
2746 : {
2747 0 : auto func = [this, &var_name](unsigned int comp) { return &coupledDot(var_name, comp); };
2748 0 : return coupledVectorHelper<const VariableValue *>(var_name, func);
2749 : }
2750 :
2751 : std::vector<const ADVariableValue *>
2752 0 : Coupleable::adCoupledDots(const std::string & var_name) const
2753 : {
2754 0 : auto func = [this, &var_name](unsigned int comp) { return &adCoupledDot(var_name, comp); };
2755 0 : return coupledVectorHelper<const ADVariableValue *>(var_name, func);
2756 : }
2757 :
2758 : template <>
2759 : const GenericVariableValue<false> &
2760 158 : Coupleable::coupledGenericDot<false>(const std::string & var_name, unsigned int comp) const
2761 : {
2762 158 : return coupledDot(var_name, comp);
2763 : }
2764 :
2765 : template <>
2766 : const GenericVariableValue<true> &
2767 56 : Coupleable::coupledGenericDot<true>(const std::string & var_name, unsigned int comp) const
2768 : {
2769 56 : return adCoupledDot(var_name, comp);
2770 : }
2771 :
2772 : // Explicit instantiations
2773 :
2774 : template const Real & Coupleable::getDefaultNodalValue<Real>(const std::string & var_name,
2775 : unsigned int comp) const;
2776 :
2777 : template const Real & Coupleable::coupledNodalValue<Real>(const std::string & var_name,
2778 : unsigned int comp) const;
2779 : template const ADReal & Coupleable::adCoupledNodalValue<Real>(const std::string & var_name,
2780 : unsigned int comp) const;
2781 : template const ADRealVectorValue &
2782 : Coupleable::adCoupledNodalValue<RealVectorValue>(const std::string & var_name,
2783 : unsigned int comp) const;
2784 :
2785 : template const RealVectorValue &
2786 : Coupleable::coupledNodalValue<RealVectorValue>(const std::string & var_name,
2787 : unsigned int comp) const;
2788 : template const Real & Coupleable::coupledNodalValueOld<Real>(const std::string & var_name,
2789 : unsigned int comp) const;
2790 : template const RealVectorValue &
2791 : Coupleable::coupledNodalValueOld<RealVectorValue>(const std::string & var_name,
2792 : unsigned int comp) const;
2793 : template const Real & Coupleable::coupledNodalValueOlder<Real>(const std::string & var_name,
2794 : unsigned int comp) const;
2795 : template const RealVectorValue &
2796 : Coupleable::coupledNodalValueOlder<RealVectorValue>(const std::string & var_name,
2797 : unsigned int comp) const;
2798 : template const Real & Coupleable::coupledNodalValuePreviousNL<Real>(const std::string & var_name,
2799 : unsigned int comp) const;
2800 : template const RealVectorValue &
2801 : Coupleable::coupledNodalValuePreviousNL<RealVectorValue>(const std::string & var_name,
2802 : unsigned int comp) const;
2803 : template const Real & Coupleable::coupledNodalDot<Real>(const std::string & var_name,
2804 : unsigned int comp) const;
2805 : template const RealVectorValue &
2806 : Coupleable::coupledNodalDot<RealVectorValue>(const std::string & var_name, unsigned int comp) const;
|