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