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