https://mooseframework.inl.gov
MooseVariableData.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "MooseVariableData.h"
11 #include "MooseVariableField.h"
12 #include "Assembly.h"
13 #include "MooseError.h"
14 #include "DisplacedSystem.h"
15 #include "TimeIntegrator.h"
16 #include "MooseVariableFE.h"
17 #include "MooseTypes.h"
18 
19 #include "libmesh/quadrature.h"
20 #include "libmesh/fe_base.h"
21 #include "libmesh/system.h"
22 #include "libmesh/type_n_tensor.h"
23 #include "libmesh/fe_interface.h"
24 
25 using namespace libMesh;
26 
27 template <typename OutputType>
29  SystemBase & sys,
30  THREAD_ID tid,
31  Moose::ElementType element_type,
32  const QBase * const & qrule_in,
33  const QBase * const & qrule_face_in,
34  const Node * const & node,
35  const Elem * const & elem)
36 
37  : MooseVariableDataBase<OutputType>(var, sys, tid),
38  _fe_type(var.feType()),
39  _var_num(var.number()),
40  _assembly(_subproblem.assembly(_tid, var.kind() == Moose::VAR_SOLVER ? sys.number() : 0)),
41  _element_type(element_type),
42  _ad_zero(0),
43  _need_ad_u_dot(false),
44  _need_ad_u_dotdot(false),
45  _need_second(false),
46  _need_second_old(false),
47  _need_second_older(false),
48  _need_second_previous_nl(false),
49  _need_curl(false),
50  _need_curl_old(false),
51  _need_curl_older(false),
52  _need_div(false),
53  _need_div_old(false),
54  _need_div_older(false),
55  _need_ad(false),
56  _need_ad_u(false),
57  _need_ad_grad_u(false),
58  _need_ad_grad_u_dot(false),
59  _need_ad_second_u(false),
60  _need_ad_curl_u(false),
61  _has_dof_indices(false),
62  _qrule(qrule_in),
63  _qrule_face(qrule_face_in),
64  _use_dual(var.useDual()),
65  _second_phi_assembly_method(nullptr),
66  _second_phi_face_assembly_method(nullptr),
67  _curl_phi_assembly_method(nullptr),
68  _curl_phi_face_assembly_method(nullptr),
69  _div_phi_assembly_method(nullptr),
70  _div_phi_face_assembly_method(nullptr),
71  _ad_grad_phi_assembly_method(nullptr),
72  _ad_grad_phi_face_assembly_method(nullptr),
73  _time_integrator(nullptr),
74  _node(node),
75  _elem(elem),
76  _displaced(dynamic_cast<const DisplacedSystem *>(&_sys) ? true : false),
77  _current_side(_assembly.side())
78 {
79  _continuity = FEInterface::get_continuity(_fe_type);
80 
82 
84 
85  // Initialize AD zero with zero derivatives
86  const auto old_do = ADReal::do_derivatives;
87  ADReal::do_derivatives = true;
88  _ad_zero = 0.;
89  ADReal::do_derivatives = old_do;
90 
91  switch (_element_type)
92  {
94  {
95  _phi_assembly_method = &Assembly::fePhi<OutputShape>;
96  _phi_face_assembly_method = &Assembly::fePhiFace<OutputShape>;
97  _grad_phi_assembly_method = &Assembly::feGradPhi<OutputShape>;
98  _grad_phi_face_assembly_method = &Assembly::feGradPhiFace<OutputShape>;
99  _second_phi_assembly_method = &Assembly::feSecondPhi<OutputShape>;
100  _second_phi_face_assembly_method = &Assembly::feSecondPhiFace<OutputShape>;
101  _curl_phi_assembly_method = &Assembly::feCurlPhi<OutputShape>;
102  _curl_phi_face_assembly_method = &Assembly::feCurlPhiFace<OutputShape>;
103  _div_phi_assembly_method = &Assembly::feDivPhi<OutputShape>;
104  _div_phi_face_assembly_method = &Assembly::feDivPhiFace<OutputShape>;
105  _ad_grad_phi_assembly_method = &Assembly::feADGradPhi<OutputShape>;
106  _ad_grad_phi_face_assembly_method = &Assembly::feADGradPhiFace<OutputShape>;
107 
110  break;
111  }
113  {
114  _phi_assembly_method = &Assembly::fePhiNeighbor<OutputShape>;
115  _phi_face_assembly_method = &Assembly::fePhiFaceNeighbor<OutputShape>;
116  _grad_phi_assembly_method = &Assembly::feGradPhiNeighbor<OutputShape>;
117  _grad_phi_face_assembly_method = &Assembly::feGradPhiFaceNeighbor<OutputShape>;
118  _second_phi_assembly_method = &Assembly::feSecondPhiNeighbor<OutputShape>;
119  _second_phi_face_assembly_method = &Assembly::feSecondPhiFaceNeighbor<OutputShape>;
120  _curl_phi_assembly_method = &Assembly::feCurlPhiNeighbor<OutputShape>;
121  _curl_phi_face_assembly_method = &Assembly::feCurlPhiFaceNeighbor<OutputShape>;
122  _div_phi_assembly_method = &Assembly::feDivPhiNeighbor<OutputShape>;
123  _div_phi_face_assembly_method = &Assembly::feDivPhiFaceNeighbor<OutputShape>;
124 
125  _ad_grad_phi = nullptr;
126  _ad_grad_phi_face = nullptr;
127  break;
128  }
130  {
131  if (_use_dual)
132  {
133  _phi_assembly_method = &Assembly::feDualPhiLower<OutputType>;
134  _phi_face_assembly_method = &Assembly::feDualPhiLower<OutputType>; // Place holder
135  _grad_phi_assembly_method = &Assembly::feGradDualPhiLower<OutputType>;
136  _grad_phi_face_assembly_method = &Assembly::feGradDualPhiLower<OutputType>; // Place holder
137  }
138  else
139  {
140  _phi_assembly_method = &Assembly::fePhiLower<OutputType>;
141  _phi_face_assembly_method = &Assembly::fePhiLower<OutputType>; // Place holder
142  _grad_phi_assembly_method = &Assembly::feGradPhiLower<OutputType>;
143  _grad_phi_face_assembly_method = &Assembly::feGradPhiLower<OutputType>; // Place holder
144  }
145 
146  _ad_grad_phi = nullptr;
147  _ad_grad_phi_face = nullptr;
148  break;
149  }
150  }
155 }
156 
157 template <typename OutputType>
158 void
160 {
161  switch (gm_type)
162  {
163  case Moose::Volume:
164  {
165  _current_qrule = _qrule;
166  _current_phi = _phi;
167  _current_grad_phi = _grad_phi;
168  _current_second_phi = _second_phi;
169  _current_curl_phi = _curl_phi;
170  _current_div_phi = _div_phi;
171  _current_ad_grad_phi = _ad_grad_phi;
172  break;
173  }
174  case Moose::Face:
175  {
176  _current_qrule = _qrule_face;
177  _current_phi = _phi_face;
178  _current_grad_phi = _grad_phi_face;
179  _current_second_phi = _second_phi_face;
180  _current_curl_phi = _curl_phi_face;
181  _current_div_phi = _div_phi_face;
182  _current_ad_grad_phi = _ad_grad_phi_face;
183  break;
184  }
185  }
186 }
187 
188 template <typename OutputType>
191 {
192  if (_sys.solutionUDot())
193  {
194  _need_u_dot = true;
195  return _u_dot;
196  }
197  else
198  mooseError("MooseVariableFE: Time derivative of solution (`u_dot`) is not stored. Please set "
199  "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
200 }
201 
202 template <typename OutputType>
205 {
206  if (_sys.solutionUDotDot())
207  {
208  _need_u_dotdot = true;
209  return _u_dotdot;
210  }
211  else
212  mooseError("MooseVariableFE: Second time derivative of solution (`u_dotdot`) is not stored. "
213  "Please set uDotDotRequested() to true in FEProblemBase before requesting "
214  "`u_dotdot`.");
215 }
216 
217 template <typename OutputType>
220 {
221  if (_sys.solutionUDotOld())
222  {
223  _need_u_dot_old = true;
224  return _u_dot_old;
225  }
226  else
227  mooseError("MooseVariableFE: Old time derivative of solution (`u_dot_old`) is not stored. "
228  "Please set uDotOldRequested() to true in FEProblemBase before requesting "
229  "`u_dot_old`.");
230 }
231 
232 template <typename OutputType>
235 {
236  if (_sys.solutionUDotDotOld())
237  {
238  _need_u_dotdot_old = true;
239  return _u_dotdot_old;
240  }
241  else
242  mooseError("MooseVariableFE: Old second time derivative of solution (`u_dotdot_old`) is not "
243  "stored. Please set uDotDotOldRequested() to true in FEProblemBase before "
244  "requesting `u_dotdot_old`");
245 }
246 
247 template <typename OutputType>
250 {
251  if (_sys.solutionUDot())
252  {
253  _need_grad_dot = true;
254  return _grad_u_dot;
255  }
256  else
257  mooseError("MooseVariableFE: Time derivative of solution (`u_dot`) is not stored. Please set "
258  "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
259 }
260 
261 template <typename OutputType>
264 {
265  if (_sys.solutionUDotDot())
266  {
267  _need_grad_dotdot = true;
268  return _grad_u_dotdot;
269  }
270  else
271  mooseError("MooseVariableFE: Second time derivative of solution (`u_dotdot`) is not stored. "
272  "Please set uDotDotRequested() to true in FEProblemBase before requesting "
273  "`u_dotdot`.");
274 }
275 
276 template <typename OutputType>
279 {
280  secondPhi();
281  secondPhiFace();
282  switch (state)
283  {
284  case Moose::Current:
285  {
286  _need_second = true;
287  return _second_u;
288  }
289 
290  case Moose::Old:
291  {
292  _need_second_old = true;
293  return _second_u_old;
294  }
295 
296  case Moose::Older:
297  {
298  _need_second_older = true;
299  return _second_u_older;
300  }
301 
302  case Moose::PreviousNL:
303  {
304  _need_second_previous_nl = true;
305  return _second_u_previous_nl;
306  }
307 
308  default:
309  // We should never get here but gcc requires that we have a default. See
310  // htpps://stackoverflow.com/questions/18680378/after-defining-case-for-all-enum-values-compiler-still-says-control-reaches-e
311  mooseError("Unknown SolutionState!");
312  }
313 }
314 
315 template <typename OutputType>
318 {
319  curlPhi();
320  curlPhiFace();
321  switch (state)
322  {
323  case Moose::Current:
324  {
325  _need_curl = true;
326  return _curl_u;
327  }
328 
329  case Moose::Old:
330  {
331  _need_curl_old = true;
332  return _curl_u_old;
333  }
334 
335  case Moose::Older:
336  {
337  _need_curl_older = true;
338  return _curl_u_older;
339  }
340 
341  default:
342  mooseError("We don't currently support curl from the previous non-linear iteration");
343  }
344 }
345 
346 template <typename OutputType>
349 {
350  divPhi();
351  divPhiFace();
352  switch (state)
353  {
354  case Moose::Current:
355  {
356  _need_div = true;
357  return _div_u;
358  }
359 
360  case Moose::Old:
361  {
362  _need_div_old = true;
363  return _div_u_old;
364  }
365 
366  case Moose::Older:
367  {
368  _need_div_older = true;
369  return _div_u_older;
370  }
371 
372  default:
373  mooseError("We don't currently support divergence from the previous non-linear iteration");
374  }
375 }
376 
377 template <typename OutputType>
380 {
381  _second_phi = &_second_phi_assembly_method(_assembly, _fe_type);
382  return *_second_phi;
383 }
384 
385 template <typename OutputType>
388 {
389  _second_phi_face = &_second_phi_face_assembly_method(_assembly, _fe_type);
390  return *_second_phi_face;
391 }
392 
393 template <typename OutputType>
396 {
397  _curl_phi = &_curl_phi_assembly_method(_assembly, _fe_type);
398  return *_curl_phi;
399 }
400 
401 template <typename OutputType>
404 {
405  _curl_phi_face = &_curl_phi_face_assembly_method(_assembly, _fe_type);
406  return *_curl_phi_face;
407 }
408 
409 template <typename OutputType>
412 {
413  _div_phi = &_div_phi_assembly_method(_assembly, _fe_type);
414  return *_div_phi;
415 }
416 
417 template <typename OutputType>
420 {
421  _div_phi_face = &_div_phi_face_assembly_method(_assembly, _fe_type);
422  return *_div_phi_face;
423 }
424 
425 template <typename OutputType>
426 template <bool monomial>
427 void
429 {
430  const auto num_dofs = _dof_indices.size();
431  if (num_dofs > 0)
432  fetchDoFValues();
433 
434  const bool is_transient = _subproblem.isTransient();
435  const auto nqp = _current_qrule->n_points();
436  auto && active_coupleable_matrix_tags = _subproblem.getActiveFEVariableCoupleableMatrixTags(_tid);
437 
438  // Map grad_phi using Eigen so that we can perform array operations easier
439  if constexpr (std::is_same_v<OutputType, RealEigenVector>)
440  {
441  if (_qrule == _current_qrule)
442  {
443  _mapped_grad_phi.resize(num_dofs);
444  for (const auto i : make_range(num_dofs))
445  {
446  _mapped_grad_phi[i].resize(nqp, Eigen::Map<RealDIMValue>(nullptr));
447  for (const auto qp : make_range(nqp))
448  // Note: this does NOT do any allocation. It is "reconstructing" the object in place
449  new (&_mapped_grad_phi[i][qp])
450  Eigen::Map<RealDIMValue>(const_cast<Real *>(&(*_current_grad_phi)[i][qp](0)));
451  }
452  }
453  else
454  {
455  _mapped_grad_phi_face.resize(num_dofs);
456  for (const auto i : make_range(num_dofs))
457  {
458  _mapped_grad_phi_face[i].resize(nqp, Eigen::Map<RealDIMValue>(nullptr));
459  for (const auto qp : make_range(nqp))
460  // Note: this does NOT do any allocation. It is "reconstructing" the object in place
461  new (&_mapped_grad_phi_face[i][qp])
462  Eigen::Map<RealDIMValue>(const_cast<Real *>(&(*_current_grad_phi)[i][qp](0)));
463  }
464  }
465  }
466 
467  mooseAssert(
468  !(_need_second || _need_second_old || _need_second_older || _need_second_previous_nl) ||
469  _current_second_phi,
470  "We're requiring a second calculation but have not set a second shape function!");
471  mooseAssert(!(_need_curl || _need_curl_old) || _current_curl_phi,
472  "We're requiring a curl calculation but have not set a curl shape function!");
473  mooseAssert(!(_need_div || _need_div_old) || _current_div_phi,
474  "We're requiring a divergence calculation but have not set a div shape function!");
475 
476  // Helper for filling values
477  const auto fill = [this, nqp, num_dofs](auto & dest, const auto & phi, const auto & dof_values)
478  {
479  if constexpr (monomial)
480  libmesh_ignore(num_dofs);
481 
482  // Deduce OutputType
483  constexpr bool is_real = std::is_same_v<OutputType, Real>;
484  constexpr bool is_real_vector = std::is_same_v<OutputType, RealVectorValue>;
485  constexpr bool is_eigen = std::is_same_v<OutputType, RealEigenVector>;
486  static_assert(is_real || is_real_vector || is_eigen, "Unsupported type");
487 
488  // this is only used in the RealEigenVector case to get this->_count
489  if constexpr (!is_eigen)
490  libmesh_ignore(this);
491 
492  // Deduce type of value within dest MooseArray
493  using dest_array_type = typename std::remove_reference_t<decltype(dest)>::value_type;
494  constexpr bool is_value = std::is_same_v<dest_array_type, OutputType>;
495  constexpr bool is_gradient = std::is_same_v<dest_array_type, OutputGradient>;
496  constexpr bool is_second = std::is_same_v<dest_array_type, OutputSecond>;
497  constexpr bool is_divergence = std::is_same_v<dest_array_type, OutputDivergence>;
498  static_assert(is_value || is_gradient || is_second || is_divergence,
499  "Unsupported destination array type");
500 
501  // Sets a value to zero at a quadrature point
502  const auto set_zero = [this, &dest](const auto qp)
503  {
504  if constexpr (!is_eigen)
505  libmesh_ignore(this);
506 
507  if constexpr (is_real || is_real_vector)
508  dest[qp] = 0;
509  else if constexpr (is_eigen)
510  {
511  if constexpr (is_value)
512  dest[qp].setZero(this->_count);
513  else if constexpr (is_gradient)
514  dest[qp].setZero(this->_count, LIBMESH_DIM);
515  else if constexpr (is_second)
516  dest[qp].setZero(this->_count, LIBMESH_DIM * LIBMESH_DIM);
517  else
518  static_assert(Moose::always_false<OutputType, dest_array_type>, "Unsupported type");
519  }
520  else
521  static_assert(Moose::always_false<OutputType, dest_array_type>, "Unsupported type");
522  };
523 
524  // Accumulates a value
525  const auto accumulate = [&dest, &phi, &dof_values](const auto i, const auto qp)
526  {
527  if constexpr (is_real || is_real_vector || (is_eigen && is_value))
528  {
529  if constexpr (is_value || is_divergence)
530  dest[qp] += phi[i][qp] * dof_values[i];
531  else if constexpr (is_gradient || is_second)
532  dest[qp].add_scaled(phi[i][qp], dof_values[i]);
533  else
534  static_assert(Moose::always_false<OutputType, dest_array_type>, "Unsupported type");
535  }
536  else if constexpr (is_eigen)
537  {
538  if constexpr (is_gradient)
539  {
540  for (const auto d : make_range(Moose::dim))
541  dest[qp].col(d) += phi[i][qp](d) * dof_values[i];
542  }
543  else if constexpr (is_second)
544  {
545  for (unsigned int d = 0, d1 = 0; d1 < LIBMESH_DIM; ++d1)
546  for (const auto d2 : make_range(Moose::dim))
547  dest[qp].col(d++) += phi[i][qp](d1, d2) * dof_values[i];
548  }
549  else
550  static_assert(Moose::always_false<OutputType, dest_array_type>, "Unsupported type");
551  }
552  else
553  static_assert(Moose::always_false<OutputType, dest_array_type>, "Unsupported type");
554  };
555 
556  dest.resize(nqp);
557 
558  // Monomial case, accumulate dest[0] and set dest[>0] to dest[0]
559  if constexpr (monomial)
560  {
561  mooseAssert(num_dofs == 1, "Should have only one dof");
562  set_zero(0);
563  accumulate(0, 0);
564  for (unsigned int qp = 1; qp < nqp; ++qp)
565  dest[qp] = dest[0];
566  }
567  // Non monomial case
568  else
569  {
570  for (const auto qp : make_range(nqp))
571  set_zero(qp);
572  for (const auto i : make_range(num_dofs))
573  for (const auto qp : make_range(nqp))
574  accumulate(i, qp);
575  }
576  };
577 
578  // Curl
579  if (_need_curl)
580  fill(_curl_u, *_current_curl_phi, _vector_tags_dof_u[_solution_tag]);
581  if (is_transient && _need_curl_old)
582  fill(_curl_u_old, *_current_curl_phi, _vector_tags_dof_u[_old_solution_tag]);
583 
584  // Div
585  if (_need_div)
586  fill(_div_u, *_current_div_phi, _vector_tags_dof_u[_solution_tag]);
587  if (is_transient && _need_div_old)
588  fill(_div_u_old, *_current_div_phi, _vector_tags_dof_u[_old_solution_tag]);
589 
590  // Second
591  if (_need_second)
592  fill(_second_u, *_current_second_phi, _vector_tags_dof_u[_solution_tag]);
593  if (_need_second_previous_nl)
594  fill(
595  _second_u_previous_nl, *_current_second_phi, _vector_tags_dof_u[_previous_nl_solution_tag]);
596 
597  // Vector tags
598  for (auto tag : _required_vector_tags)
599  {
600  if (_need_vector_tag_u[tag] && _sys.hasVector(tag))
601  {
602  mooseAssert(_sys.getVector(tag).closed(), "Vector should be closed");
603  fill(_vector_tag_u[tag], *_current_phi, _vector_tags_dof_u[tag]);
604  }
605  if (_need_vector_tag_grad[tag] && _sys.hasVector(tag))
606  {
607  mooseAssert(_sys.getVector(tag).closed(), "Vector should be closed");
608  fill(_vector_tag_grad[tag], *_current_grad_phi, _vector_tags_dof_u[tag]);
609  }
610  }
611 
612  // Matrix tags
613  for (auto tag : active_coupleable_matrix_tags)
614  if (_need_matrix_tag_u[tag])
615  fill(_matrix_tag_u[tag], *_current_phi, _matrix_tags_dof_u[tag]);
616 
617  // Derivatives and old values
618  if (is_transient)
619  {
620  if (_need_second_old)
621  fill(_second_u_old, *_current_second_phi, _vector_tags_dof_u[_old_solution_tag]);
622  if (_need_second_older)
623  fill(_second_u_older, *_current_second_phi, _vector_tags_dof_u[_older_solution_tag]);
624  if (_need_u_dot)
625  fill(_u_dot, *_current_phi, _dof_values_dot);
626  if (_need_u_dotdot)
627  fill(_u_dotdot, *_current_phi, _dof_values_dotdot);
628  if (_need_u_dot_old)
629  fill(_u_dot_old, *_current_phi, _dof_values_dot_old);
630  if (_need_u_dotdot_old)
631  fill(_u_dotdot_old, *_current_phi, _dof_values_dotdot_old);
632 
633  if (_need_du_dot_du)
634  {
635  _du_dot_du.resize(nqp);
636  for (const auto i : make_range(num_dofs))
637  for (const auto qp : make_range(nqp))
638  _du_dot_du[qp] = _dof_du_dot_du[i];
639  }
640  if (_need_du_dotdot_du)
641  {
642  _du_dotdot_du.resize(nqp);
643  for (const auto i : make_range(num_dofs))
644  for (const auto qp : make_range(nqp))
645  _du_dotdot_du[qp] = _dof_du_dotdot_du[i];
646  }
647 
648  if (_need_grad_dot)
649  fill(_grad_u_dot, *_current_grad_phi, _dof_values_dot);
650  if (_need_grad_dotdot)
651  fill(_grad_u_dotdot, *_current_grad_phi, _dof_values_dotdot);
652  }
653 
654  // Automatic differentiation, not for eigen
655  if constexpr (!std::is_same_v<OutputType, RealEigenVector>)
656  if (_need_ad)
657  computeAD(num_dofs, nqp);
658 }
659 
660 template <typename OutputType>
661 void
663 {
664  computeValuesInternal</* monomial = */ false>();
665 }
666 
667 template <typename OutputType>
668 void
670 {
671  if (_dof_indices.size() == 0)
672  return;
673 
674  // Monomial optimizations are not appropriate after p-refinement
675  if (_elem->p_level())
676  computeValues();
677  else
678  computeValuesInternal</* monomial = */ true>();
679 }
680 
681 template <typename OutputType>
682 void
683 MooseVariableData<OutputType>::computeAD(const unsigned int num_dofs, const unsigned int nqp)
684 {
685  const bool do_derivatives = Moose::doDerivatives(_subproblem, _sys);
686 
687  _ad_dof_values.resize(num_dofs);
688  for (const auto i : make_range(num_dofs))
689  _ad_dof_values[i] = _vector_tags_dof_u[_solution_tag][i];
690  // NOTE! You have to do this AFTER setting the value!
691  if (do_derivatives)
692  for (const auto i : make_range(num_dofs))
693  Moose::derivInsert(_ad_dof_values[i].derivatives(), _dof_indices[i], 1.);
694 
695  if (_need_ad_u)
696  {
697  _ad_u.resize(nqp);
698  for (const auto qp : make_range(nqp))
699  _ad_u[qp] = _ad_zero;
700 
701  for (const auto i : make_range(num_dofs))
702  for (const auto qp : make_range(nqp))
703  _ad_u[qp] += _ad_dof_values[i] * (*_current_phi)[i][qp];
704  }
705 
706  if (_need_ad_grad_u)
707  {
708  _ad_grad_u.resize(nqp);
709  for (const auto qp : make_range(nqp))
710  _ad_grad_u[qp] = _ad_zero;
711 
712  // The latter check here is for handling the fact that we have not yet implemented
713  // calculation of ad_grad_phi for neighbor and neighbor-face, so if we are in that
714  // situation we need to default to using the non-ad grad_phi
715  if (_displaced && _current_ad_grad_phi)
716  for (const auto i : make_range(num_dofs))
717  for (const auto qp : make_range(nqp))
718  _ad_grad_u[qp] += _ad_dof_values[i] * (*_current_ad_grad_phi)[i][qp];
719  else
720  for (const auto i : make_range(num_dofs))
721  for (const auto qp : make_range(nqp))
722  _ad_grad_u[qp] += _ad_dof_values[i] * (*_current_grad_phi)[i][qp];
723  }
724 
725  if (_need_ad_second_u)
726  {
727  _ad_second_u.resize(nqp);
728  for (const auto qp : make_range(nqp))
729  _ad_second_u[qp] = _ad_zero;
730 
731  for (const auto i : make_range(num_dofs))
732  for (const auto qp : make_range(nqp))
733  // Note that this will not carry any derivatives with respect to displacements because
734  // those calculations have not yet been implemented in Assembly
735  _ad_second_u[qp] += _ad_dof_values[i] * (*_current_second_phi)[i][qp];
736  }
737 
738  if (_need_ad_curl_u)
739  {
740  _ad_curl_u.resize(nqp);
741  for (const auto qp : make_range(nqp))
742  _ad_curl_u[qp] = _ad_zero;
743 
744  for (const auto i : make_range(num_dofs))
745  for (const auto qp : make_range(nqp))
746  {
747  mooseAssert(_current_curl_phi,
748  "We're requiring a curl calculation but have not set a curl shape function!");
749 
750  // Note that the current version of _ad_curl_u is not yet implemented for mesh displacement
751  _ad_curl_u[qp] += _ad_dof_values[i] * (*_current_curl_phi)[i][qp];
752  }
753  }
754 
755  bool is_transient = _subproblem.isTransient();
756  if (is_transient)
757  {
758  if (_need_ad_u_dot)
759  {
760  _ad_dofs_dot.resize(num_dofs);
761  if (_need_ad_u_dotdot)
762  _ad_dofs_dotdot.resize(num_dofs);
763  _ad_u_dot.resize(nqp);
764 
765  if (_time_integrator && _time_integrator->dt())
766  {
767  for (const auto i : make_range(num_dofs))
768  _ad_dofs_dot[i] = _ad_dof_values[i];
769  for (const auto i : make_range(num_dofs))
770  _time_integrator->computeADTimeDerivatives(_ad_dofs_dot[i],
771  _dof_indices[i],
772  _need_ad_u_dotdot ? _ad_dofs_dotdot[i]
773  : _ad_real_dummy);
774 
775  for (const auto qp : make_range(nqp))
776  _ad_u_dot[qp] = _ad_zero;
777  for (const auto i : make_range(num_dofs))
778  for (const auto qp : make_range(nqp))
779  _ad_u_dot[qp] += (*_current_phi)[i][qp] * _ad_dofs_dot[i];
780  }
781  // We are too early in the setup to have a time integrator, so we are not really using the
782  // AD-derivatives. We set the AD value of the derivatives to the nonAD value
783  else if (!_time_integrator)
784  {
785  for (const auto i : make_range(num_dofs))
786  _ad_dofs_dot[i] = _dof_values_dot[i];
787  for (const auto qp : make_range(nqp))
788  _ad_u_dot[qp] = _u_dot[qp];
789  }
790  }
791 
792  if (_need_ad_u_dotdot)
793  {
794  _ad_u_dotdot.resize(nqp);
795  for (const auto qp : make_range(nqp))
796  _ad_u_dotdot[qp] = _ad_zero;
797 
798  if (_time_integrator && _time_integrator->dt())
799  for (const auto i : make_range(num_dofs))
800  for (const auto qp : make_range(nqp))
801  _ad_u_dotdot[qp] += (*_current_phi)[i][qp] * _ad_dofs_dotdot[i];
802  else if (!_time_integrator)
803  for (const auto qp : make_range(nqp))
804  _ad_u_dotdot[qp] = _u_dotdot[qp];
805  }
806 
807  if (_need_ad_grad_u_dot)
808  {
809  _ad_grad_u_dot.resize(nqp);
810 
811  if (_time_integrator && _time_integrator->dt())
812  {
813  for (const auto qp : make_range(nqp))
814  _ad_grad_u_dot[qp] = _ad_zero;
815 
816  // The latter check here is for handling the fact that we have not yet implemented
817  // calculation of ad_grad_phi for neighbor and neighbor-face, so if we are in that
818  // situation we need to default to using the non-ad grad_phi
819  if (_displaced && _current_ad_grad_phi)
820  for (const auto i : make_range(num_dofs))
821  for (const auto qp : make_range(nqp))
822  _ad_grad_u_dot[qp] += _ad_dofs_dot[i] * (*_current_ad_grad_phi)[i][qp];
823  else
824  for (const auto i : make_range(num_dofs))
825  for (const auto qp : make_range(nqp))
826  _ad_grad_u_dot[qp] += _ad_dofs_dot[i] * (*_current_grad_phi)[i][qp];
827  }
828  else if (!_time_integrator)
829  for (const auto qp : make_range(nqp))
830  _ad_grad_u_dot[qp] = _grad_u_dot[qp];
831  }
832  }
833 }
834 
835 template <>
836 void
837 MooseVariableData<RealEigenVector>::computeAD(const unsigned int /*num_dofs*/,
838  const unsigned int /*nqp*/)
839 {
840  mooseError("AD for array variable has not been implemented");
841 }
842 
843 template <typename OutputType>
844 void
845 MooseVariableData<OutputType>::setDofValue(const OutputData & value, unsigned int index)
846 {
847  auto & dof_values = _vector_tags_dof_u[_solution_tag];
848  dof_values[index] = value;
849  _has_dof_values = true;
850 
851  auto & u = _vector_tag_u[_solution_tag];
852  const auto nqps = u.size();
853  const auto ndofs = dof_values.size();
854  for (const auto qp : make_range(nqps))
855  u[qp] *= 0.;
856  for (const auto qp : make_range(nqps))
857  for (const auto i : make_range(ndofs))
858  u[qp] += (*_phi)[i][qp] * dof_values[i];
859 }
860 
861 template <typename OutputType>
862 void
864 {
865  auto & dof_values = _vector_tags_dof_u[_solution_tag];
866  for (unsigned int i = 0; i < values.size(); i++)
867  dof_values[i] = values(i);
868 
869  _has_dof_values = true;
870 
871  auto & u = _vector_tag_u[_solution_tag];
872  const auto nqps = u.size();
873  const auto ndofs = dof_values.size();
874  for (const auto qp : make_range(nqps))
875  u[qp] *= 0.;
876  for (const auto qp : make_range(nqps))
877  for (const auto i : make_range(ndofs))
878  u[qp] += (*_phi)[i][qp] * dof_values[i];
879 }
880 
881 template <typename OutputType>
882 void
884  const OutputData & v)
885 {
886  residual.set(_nodal_dof_index, v);
887 }
888 
889 template <>
890 void
892  const RealEigenVector & v)
893 {
894  for (const auto j : make_range(_count))
895  residual.set(_nodal_dof_index + j, v(j));
896 }
897 
898 template <typename OutputType>
901 {
902  mooseAssert(_subproblem.mesh().isSemiLocal(const_cast<Node *>(&node)), "Node is not Semilocal");
903 
904  // Make sure that the node has DOFs
905  /* Note, this is a reproduction of an assert within libMesh::Node::dof_number, this is done to
906  * produce a better error (see misc/check_error.node_value_off_block) */
907  mooseAssert(node.n_dofs(_sys.number(), _var_num) > 0,
908  "Node " << node.id() << " does not contain any dofs for the "
909  << _sys.system().variable_name(_var_num) << " variable");
910 
911  dof_id_type dof = node.dof_number(_sys.number(), _var_num, 0);
912 
913  switch (state)
914  {
915  case Moose::Current:
916  return (*_sys.currentSolution())(dof);
917 
918  case Moose::Old:
919  return _sys.solutionOld()(dof);
920 
921  case Moose::Older:
922  return _sys.solutionOlder()(dof);
923 
924  default:
925  mooseError("PreviousNL not currently supported for getNodalValue");
926  }
927 }
928 
929 template <>
932  Moose::SolutionState state) const
933 {
934  mooseAssert(_subproblem.mesh().isSemiLocal(const_cast<Node *>(&node)), "Node is not Semilocal");
935 
936  // Make sure that the node has DOFs
937  /* Note, this is a reproduction of an assert within libMesh::Node::dof_number, this is done to
938  * produce a better error (see misc/check_error.node_value_off_block) */
939  mooseAssert(node.n_dofs(_sys.number(), _var_num) > 0,
940  "Node " << node.id() << " does not contain any dofs for the "
941  << _sys.system().variable_name(_var_num) << " variable");
942 
943  dof_id_type dof = node.dof_number(_sys.number(), _var_num, 0);
944 
945  RealEigenVector v(_count);
946  switch (state)
947  {
948  case Moose::Current:
949  for (unsigned int i = 0; i < _count; ++i)
950  v(i) = (*_sys.currentSolution())(dof++);
951  break;
952 
953  case Moose::Old:
954  for (unsigned int i = 0; i < _count; ++i)
955  v(i) = _sys.solutionOld()(dof++);
956  break;
957 
958  case Moose::Older:
959  for (unsigned int i = 0; i < _count; ++i)
960  v(i) = _sys.solutionOlder()(dof++);
961  break;
962 
963  default:
964  mooseError("PreviousNL not currently supported for getNodalValue");
965  }
966  return v;
967 }
968 
969 template <typename OutputType>
972  Moose::SolutionState state,
973  unsigned int idx) const
974 {
975  std::vector<dof_id_type> dof_indices;
976  _dof_map.dof_indices(elem, dof_indices, _var_num);
977 
978  switch (state)
979  {
980  case Moose::Current:
981  return (*_sys.currentSolution())(dof_indices[idx]);
982 
983  case Moose::Old:
984  return _sys.solutionOld()(dof_indices[idx]);
985 
986  case Moose::Older:
987  return _sys.solutionOlder()(dof_indices[idx]);
988 
989  default:
990  mooseError("PreviousNL not currently supported for getElementalValue");
991  }
992 }
993 
994 template <>
997  Moose::SolutionState state,
998  unsigned int idx) const
999 {
1000  std::vector<dof_id_type> dof_indices;
1001  _dof_map.dof_indices(elem, dof_indices, _var_num);
1002 
1003  dof_id_type dof = dof_indices[idx];
1004 
1005  RealEigenVector v(_count);
1006 
1007  switch (state)
1008  {
1009  case Moose::Current:
1010  for (unsigned int i = 0; i < _count; ++i)
1011  v(i) = (*_sys.currentSolution())(dof++);
1012  break;
1013 
1014  case Moose::Old:
1015  for (unsigned int i = 0; i < _count; ++i)
1016  v(i) = _sys.solutionOld()(dof++);
1017  break;
1018 
1019  case Moose::Older:
1020  for (unsigned int i = 0; i < _count; ++i)
1021  v(i) = _sys.solutionOlder()(dof++);
1022  break;
1023 
1024  default:
1025  mooseError("PreviousNL not currently supported for getElementalValue");
1026  }
1027  return v;
1028 }
1029 
1030 template <typename OutputType>
1031 void
1033  std::vector<dof_id_type> & dof_indices) const
1034 {
1035  _dof_map.dof_indices(elem, dof_indices, _var_num);
1036 }
1037 
1038 template <typename OutputType>
1039 void
1041  const DenseVector<Number> & v) const
1042 {
1043  sol.add_vector(v, _dof_indices);
1044 }
1045 
1046 template <>
1047 void
1049  const DenseVector<Number> & v) const
1050 {
1051  unsigned int p = 0;
1052  const auto num_dofs = _dof_indices.size();
1053  for (const auto j : make_range(_count))
1054  {
1055  unsigned int inc = (isNodal() ? j : j * num_dofs);
1056  for (const auto i : make_range(num_dofs))
1057  sol.add(_dof_indices[i] + inc, v(p++));
1058  }
1059 }
1060 
1061 template <typename OutputType>
1064 {
1065  if (_sys.solutionUDot())
1066  {
1067  _need_dof_values_dot = true;
1068  return _dof_values_dot;
1069  }
1070  else
1071  mooseError("MooseVariableData: Time derivative of solution (`u_dot`) is not stored. Please set "
1072  "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
1073 }
1074 
1075 template <typename OutputType>
1078 {
1079  if (_sys.solutionUDotDot())
1080  {
1081  _need_dof_values_dotdot = true;
1082  return _dof_values_dotdot;
1083  }
1084  else
1085  mooseError("MooseVariableData: Second time derivative of solution (`u_dotdot`) is not stored. "
1086  "Please set uDotDotRequested() to true in FEProblemBase before requesting "
1087  "`u_dotdot`.");
1088 }
1089 
1090 template <typename OutputType>
1093 {
1094  if (_sys.solutionUDotOld())
1095  {
1096  _need_dof_values_dot_old = true;
1097  return _dof_values_dot_old;
1098  }
1099  else
1100  mooseError("MooseVariableData: Old time derivative of solution (`u_dot_old`) is not stored. "
1101  "Please set uDotOldRequested() to true in FEProblemBase before requesting "
1102  "`u_dot_old`.");
1103 }
1104 
1105 template <typename OutputType>
1108 {
1109  if (_sys.solutionUDotDotOld())
1110  {
1111  _need_dof_values_dotdot_old = true;
1112  return _dof_values_dotdot_old;
1113  }
1114  else
1115  mooseError("MooseVariableData: Old second time derivative of solution (`u_dotdot_old`) is not "
1116  "stored. Please set uDotDotOldRequested() to true in FEProblemBase before "
1117  "requesting `u_dotdot_old`.");
1118 }
1119 
1120 template <typename OutputType>
1121 const MooseArray<Number> &
1123 {
1124  _need_dof_du_dot_du = true;
1125  return _dof_du_dot_du;
1126 }
1127 
1128 template <typename OutputType>
1129 const MooseArray<Number> &
1131 {
1132  _need_dof_du_dotdot_du = true;
1133  return _dof_du_dotdot_du;
1134 }
1135 
1136 template <typename OutputType>
1137 void
1139 {
1140  unsigned int nqp = _qrule->n_points();
1141 
1142  _increment.resize(nqp);
1143  // Compute the increment at each quadrature point
1144  unsigned int num_dofs = _dof_indices.size();
1145  for (const auto qp : make_range(nqp))
1146  {
1147  _increment[qp] = 0.;
1148  for (const auto i : make_range(num_dofs))
1149  _increment[qp] += (*_phi)[i][qp] * increment_vec(_dof_indices[i]);
1150  }
1151 }
1152 
1153 template <>
1154 void
1156  const NumericVector<Number> & increment_vec)
1157 {
1158  unsigned int nqp = _qrule->n_points();
1159 
1160  _increment.resize(nqp);
1161  // Compute the increment at each quadrature point
1162  unsigned int num_dofs = _dof_indices.size();
1163  if (isNodal())
1164  {
1165  for (const auto qp : make_range(nqp))
1166  {
1167  for (const auto i : make_range(num_dofs))
1168  for (const auto j : make_range(_count))
1169  _increment[qp](j) += (*_phi)[i][qp] * increment_vec(_dof_indices[i] + j);
1170  }
1171  }
1172  else
1173  {
1174  for (const auto qp : make_range(nqp))
1175  {
1176  unsigned int n = 0;
1177  for (const auto j : make_range(_count))
1178  for (const auto i : make_range(num_dofs))
1179  {
1180  _increment[qp](j) += (*_phi)[i][qp] * increment_vec(_dof_indices[i] + n);
1181  n += num_dofs;
1182  }
1183  }
1184  }
1185 }
1186 
1187 template <typename OutputType>
1188 void
1190 {
1191  if (!isNodal())
1192  mooseError("computeIncrementAtNode can only be called for nodal variables");
1193 
1194  _increment.resize(1);
1195 
1196  // Compute the increment for the current DOF
1197  _increment[0] = increment_vec(_dof_indices[0]);
1198 }
1199 
1200 template <>
1201 void
1203  const NumericVector<Number> & increment_vec)
1204 {
1205  if (!isNodal())
1206  mooseError("computeIncrementAtNode can only be called for nodal variables");
1207 
1208  _increment.resize(1);
1209 
1210  // Compute the increment for the current DOF
1211  if (isNodal())
1212  for (unsigned int j = 0; j < _count; j++)
1213  _increment[0](j) = increment_vec(_dof_indices[0] + j);
1214  else
1215  {
1216  unsigned int n = 0;
1217  const auto n_dof_indices = _dof_indices.size();
1218  for (const auto j : make_range(_count))
1219  {
1220  _increment[0](j) = increment_vec(_dof_indices[0] + n);
1221  n += n_dof_indices;
1222  }
1223  }
1224 }
1225 
1226 template <typename OutputType>
1227 const OutputType &
1229 {
1230  if (isNodal())
1231  {
1232  if (_sys.solutionUDot())
1233  {
1234  _need_dof_values_dot = true;
1235  return _nodal_value_dot;
1236  }
1237  else
1238  mooseError(
1239  "MooseVariableData: Time derivative of solution (`u_dot`) is not stored. Please set "
1240  "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
1241  }
1242  else
1243  mooseError("Nodal values can be requested only on nodal variables, variable '",
1244  var().name(),
1245  "' is not nodal.");
1246 }
1247 
1248 template <typename OutputType>
1249 const OutputType &
1251 {
1252  if (isNodal())
1253  {
1254  if (_sys.solutionUDotDot())
1255  {
1256  _need_dof_values_dotdot = true;
1257  return _nodal_value_dotdot;
1258  }
1259  else
1260  mooseError(
1261  "MooseVariableData: Second time derivative of solution (`u_dotdot`) is not stored. "
1262  "Please set uDotDotRequested() to true in FEProblemBase before requesting "
1263  "`u_dotdot`.");
1264  }
1265  else
1266  mooseError("Nodal values can be requested only on nodal variables, variable '",
1267  var().name(),
1268  "' is not nodal.");
1269 }
1270 
1271 template <typename OutputType>
1272 const OutputType &
1274 {
1275  if (isNodal())
1276  {
1277  if (_sys.solutionUDotOld())
1278  {
1279  _need_dof_values_dot_old = true;
1280  return _nodal_value_dot_old;
1281  }
1282  else
1283  mooseError("MooseVariableData: Old time derivative of solution (`u_dot_old`) is not stored. "
1284  "Please set uDotOldRequested() to true in FEProblemBase before requesting "
1285  "`u_dot_old`.");
1286  }
1287  else
1288  mooseError("Nodal values can be requested only on nodal variables, variable '",
1289  var().name(),
1290  "' is not nodal.");
1291 }
1292 
1293 template <typename OutputType>
1294 const OutputType &
1296 {
1297  if (isNodal())
1298  {
1299  if (_sys.solutionUDotDotOld())
1300  {
1301  _need_dof_values_dotdot_old = true;
1302  return _nodal_value_dotdot_old;
1303  }
1304  else
1305  mooseError(
1306  "MooseVariableData: Old second time derivative of solution (`u_dotdot_old`) is not "
1307  "stored. Please set uDotDotOldRequested() to true in FEProblemBase before "
1308  "requesting `u_dotdot_old`.");
1309  }
1310  else
1311  mooseError("Nodal values can be requested only on nodal variables, variable '",
1312  var().name(),
1313  "' is not nodal.");
1314 }
1315 
1316 template <typename OutputType>
1317 void
1319 {
1320  if (_has_dof_indices)
1321  {
1322  fetchDoFValues();
1323  assignNodalValue();
1324 
1325  if (_need_ad)
1326  fetchADNodalValues();
1327  }
1328  else
1329  zeroSizeDofValues();
1330 }
1331 
1332 template <typename OutputType>
1333 void
1335 {
1336  auto n = _dof_indices.size();
1337  libmesh_assert(n);
1338  _ad_dof_values.resize(n);
1339 
1340  if (_need_ad_u_dot)
1341  _ad_dofs_dot.resize(n);
1342  if (_need_ad_u_dotdot)
1343  _ad_dofs_dotdot.resize(n);
1344 
1345  const bool do_derivatives =
1346  ADReal::do_derivatives && _sys.number() == _subproblem.currentNlSysNum();
1347 
1348  for (decltype(n) i = 0; i < n; ++i)
1349  {
1350  _ad_dof_values[i] = _vector_tags_dof_u[_solution_tag][i];
1351  if (do_derivatives)
1352  Moose::derivInsert(_ad_dof_values[i].derivatives(), _dof_indices[i], 1.);
1353  assignADNodalValue(_ad_dof_values[i], i);
1354 
1355  if (_need_ad_u_dot)
1356  {
1357  if (_time_integrator && _time_integrator->dt())
1358  {
1359  _ad_dofs_dot[i] = _ad_dof_values[i];
1360  _time_integrator->computeADTimeDerivatives(_ad_dofs_dot[i],
1361  _dof_indices[i],
1362  _need_ad_u_dotdot ? _ad_dofs_dotdot[i]
1363  : _ad_real_dummy);
1364  }
1365  // Executing something with a time derivative at initial should not put a NaN
1366  else if (_time_integrator && !_time_integrator->dt())
1367  _ad_dofs_dot[i] = 0.;
1368  else
1369  mooseError("AD nodal time derivatives not implemented for variables without a time "
1370  "integrator (auxiliary variables)");
1371  }
1372  }
1373 }
1374 
1375 template <>
1376 void
1378 {
1379  mooseError("I do not know how to support AD with array variables");
1380 }
1381 
1382 template <>
1383 void
1384 MooseVariableData<Real>::assignADNodalValue(const ADReal & value, const unsigned int &)
1385 {
1386  _ad_nodal_value = value;
1387 }
1388 
1389 template <>
1390 void
1392  const unsigned int & component)
1393 {
1394  _ad_nodal_value(component) = value;
1395 }
1396 
1397 template <typename OutputType>
1398 void
1400 {
1401  _dof_map.dof_indices(_elem, _dof_indices, _var_num);
1402  _vector_tags_dof_u[_solution_tag].resize(_dof_indices.size());
1403 
1404  unsigned int nqp = _qrule->n_points();
1405  _vector_tag_u[_solution_tag].resize(nqp);
1406 }
1407 
1408 template <typename OutputType>
1409 void
1411 {
1412  _dof_map.dof_indices(_elem, _dof_indices, _var_num);
1413  _has_dof_values = false;
1414 
1415  // FIXME: remove this when the Richard's module is migrated to use the new NodalCoupleable
1416  // interface.
1417  _has_dof_indices = _dof_indices.size();
1418 }
1419 
1420 template <typename OutputType>
1421 void
1423 {
1424  if (std::size_t n_dofs = _node->n_dofs(_sys.number(), _var_num))
1425  {
1426  _dof_indices.resize(n_dofs);
1427  for (std::size_t i = 0; i < n_dofs; ++i)
1428  _dof_indices[i] = _node->dof_number(_sys.number(), _var_num, i);
1429  // For standard variables. _nodal_dof_index is retrieved by nodalDofIndex() which is used in
1430  // NodalBC for example
1431  _nodal_dof_index = _dof_indices[0];
1432  _has_dof_indices = true;
1433  }
1434  else
1435  {
1436  _dof_indices.clear(); // Clear these so Assembly doesn't think there's dofs here
1437  _has_dof_indices = false;
1438  }
1439 }
1440 
1441 template <typename OutputType>
1442 void
1444 {
1445  /* FIXME: this method is only for elemental auxiliary variables, so
1446  * we may want to rename it */
1447  if (_elem)
1448  {
1449  _dof_map.dof_indices(_elem, _dof_indices, _var_num);
1450  if (_elem->n_dofs(_sys.number(), _var_num) > 0)
1451  {
1452  _nodal_dof_index = _dof_indices[0];
1453 
1454  fetchDoFValues();
1455 
1456  const auto num_dofs = _dof_indices.size();
1457  for (auto & dof_u : _vector_tags_dof_u)
1458  dof_u.resize(num_dofs);
1459 
1460  for (auto & dof_u : _matrix_tags_dof_u)
1461  dof_u.resize(num_dofs);
1462 
1463  _has_dof_indices = true;
1464  }
1465  else
1466  _has_dof_indices = false;
1467  }
1468  else
1469  _has_dof_indices = false;
1470 }
1471 
1472 template <typename OutputType>
1473 void
1474 MooseVariableData<OutputType>::reinitNodes(const std::vector<dof_id_type> & nodes)
1475 {
1476  _dof_indices.clear();
1477  for (const auto & node_id : nodes)
1478  {
1479  auto && nd = _subproblem.mesh().getMesh().query_node_ptr(node_id);
1480  if (nd && (_subproblem.mesh().isSemiLocal(const_cast<Node *>(nd))))
1481  {
1482  if (nd->n_dofs(_sys.number(), _var_num) > 0)
1483  {
1484  dof_id_type dof = nd->dof_number(_sys.number(), _var_num, 0);
1485  _dof_indices.push_back(dof);
1486  }
1487  }
1488  }
1489 
1490  if (!_dof_indices.empty())
1491  _has_dof_indices = true;
1492  else
1493  _has_dof_indices = false;
1494 }
1495 
1496 template class MooseVariableData<Real>;
1497 template class MooseVariableData<RealVectorValue>;
1498 template class MooseVariableData<RealEigenVector>;
const FieldVariableDivergence & divSln(Moose::SolutionState state) const
Local solution divergence getter.
std::string name(const ElemQuality q)
const Assembly & _assembly
OutputData getNodalValue(const Node &node, Moose::SolutionState state) const
const FieldVariableValue & uDotDotOld() const
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
void computeAD(const unsigned int num_dofs, const unsigned int nqp)
compute AD things
const FieldVariablePhiCurl & curlPhi() const
curl_phi getter
void reinitNodes(const std::vector< dof_id_type > &nodes)
Set _dof_indices to the degrees of freedom existing on the passed-in nodes.
void prepare()
Get the dof indices corresponding to the current element.
Class for stuff related to variables.
void setDofValue(const OutputData &value, unsigned int index)
dof value setters
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
void computeValues()
compute the variable values
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
const FieldVariableValue & uDotDot() const
void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const
std::function< const typename OutputTools< OutputShape >::VariablePhiCurl &(const Assembly &, libMesh::FEType)> _curl_phi_assembly_method
Moose::DOFType< OutputType >::type OutputData
const OutputType & nodalValueDotOld() const
std::function< const typename OutputTools< OutputShape >::VariablePhiCurl &(const Assembly &, libMesh::FEType)> _curl_phi_face_assembly_method
void computeIncrementAtNode(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change at the current node based on increment_vec.
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:162
const FieldVariableValue & uDot() const
void computeMonomialValues()
compute the values for const monomial variables
MooseVariableData(const MooseVariableField< OutputType > &var, SystemBase &sys, THREAD_ID tid, Moose::ElementType element_type, const QBase *const &qrule_in, const QBase *const &qrule_face_in, const Node *const &node, const Elem *const &elem)
void computeNodalValues()
compute nodal things
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
ElementType
Definition: MooseTypes.h:764
std::function< const ADTemplateVariablePhiGradient< OutputShape > &(const Assembly &, libMesh::FEType)> _ad_grad_phi_assembly_method
Base class for a system (of equations)
Definition: SystemBase.h:84
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
const unsigned int _var_num
void prepareIC()
prepare the initial condition
void assignADNodalValue(const ADReal &value, const unsigned int &component)
Helper methods for assigning nodal values from their corresponding solution values (dof values as the...
const FieldVariableGradient & gradSlnDot() const
Local time derivative of solution gradient getter.
bool doDerivatives(const SubProblem &subproblem, const SystemBase &sys)
Definition: ADUtils.C:83
void setDofValues(const DenseVector< OutputData > &values)
Set local DOF values and evaluate the values on quadrature points.
unsigned int n_dofs(const unsigned int s, const unsigned int var=libMesh::invalid_uint) const
std::function< const typename OutputTools< OutputShape >::VariablePhiValue &(const Assembly &, libMesh::FEType)> _phi_face_assembly_method
OutputData getElementalValue(const Elem *elem, Moose::SolutionState state, unsigned int idx=0) const
void insertNodalValue(libMesh::NumericVector< libMesh::Number > &residual, const OutputData &v)
Write a nodal value to the passed-in solution vector.
std::function< const typename OutputTools< OutputShape >::VariablePhiSecond &(const Assembly &, libMesh::FEType)> _second_phi_face_assembly_method
void libmesh_ignore(const Args &...)
const DoFValue & dofValuesDotDotOld() const
const MooseArray< libMesh::Number > & dofValuesDuDotDu() const
const DoFValue & dofValuesDot() const
const FieldVariablePhiDivergence & divPhi() const
divergence_phi getter
dof_id_type id() const
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
void computeValuesInternal()
Internal method for computeValues() and computeMonomialValues()
const OutputType & nodalValueDot() const
const FieldVariablePhiSecond & secondPhiFace() const
second_phi_face getter
void computeIncrementAtQps(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change in solution at QPs based on increment_vec.
C_ZERO
const FieldVariableGradient & gradSlnDotDot() const
Local second time derivative of solution gradient getter.
std::function< const typename OutputTools< OutputShape >::VariablePhiDivergence &(const Assembly &, libMesh::FEType)> _div_phi_face_assembly_method
libmesh_assert(ctx)
void reinitNode()
Prepare degrees of freedom for the current node.
std::function< const ADTemplateVariablePhiGradient< OutputShape > &(const Assembly &, libMesh::FEType)> _ad_grad_phi_face_assembly_method
const TimeIntegrator * _time_integrator
Pointer to time integrator.
const DoFValue & dofValuesDotDot() const
libMesh::FEContinuity _continuity
Continuity type of the variable.
Moose::ElementType _element_type
The element type this object is storing data for. This is either Element, Neighbor, or Lower.
ADReal _ad_zero
A zero AD variable.
std::function< const typename OutputTools< OutputShape >::VariablePhiGradient &(const Assembly &, libMesh::FEType)> _grad_phi_assembly_method
const ADTemplateVariablePhiGradient< OutputShape > * _ad_grad_phi
const OutputType & nodalValueDotDotOld() const
std::function< const typename OutputTools< OutputShape >::VariablePhiGradient &(const Assembly &, libMesh::FEType)> _grad_phi_face_assembly_method
const FieldVariablePhiValue * _phi_face
std::function< const typename OutputTools< OutputShape >::VariablePhiSecond &(const Assembly &, libMesh::FEType)> _second_phi_assembly_method
GeometryType
Definition: MooseTypes.h:249
const FieldVariableCurl & curlSln(Moose::SolutionState state) const
Local solution curl getter.
C_ONE
const FieldVariableValue & uDotOld() const
SolutionState
Definition: MooseTypes.h:233
const libMesh::FEType & _fe_type
const ADTemplateVariablePhiGradient< OutputShape > * _ad_grad_phi_face
const FieldVariablePhiSecond & secondPhi() const
second_phi getter
const OutputType & nodalValueDotDot() const
SystemBase & _sys
The MOOSE system which ultimately holds the vectors and matrices relevant to this variable data...
IntRange< T > make_range(T beg, T end)
virtual unsigned int size() const override final
const TimeIntegrator * queryTimeIntegrator(const unsigned int var_num) const
Retrieve the time integrator that integrates the given variable&#39;s equation.
Definition: SystemBase.C:1673
const FieldVariablePhiGradient * _grad_phi_face
void addSolution(libMesh::NumericVector< libMesh::Number > &sol, const DenseVector< libMesh::Number > &v) const
Add passed in local DOF values to a solution vector.
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
Definition: ADReal.h:21
const FieldVariablePhiCurl & curlPhiFace() const
curl_phi_face getter
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
virtual void set(const numeric_index_type i, const T value)=0
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
const DoFValue & dofValuesDotOld() const
void setGeometry(Moose::GeometryType gm_type)
Set the geometry type before calculating variables values.
void reinitAux()
Prepare dof indices and solution values for elemental auxiliary variables.
virtual void add(const numeric_index_type i, const T value)=0
std::function< const typename OutputTools< OutputType >::VariablePhiValue &(const Assembly &, libMesh::FEType)> _phi_assembly_method
const FieldVariablePhiDivergence & divPhiFace() const
divergence_phi_face getter
const FieldVariableSecond & secondSln(Moose::SolutionState state) const
Local solution second spatial derivative getter.
bool _is_nodal
if variable is nodal
const FieldVariablePhiGradient * _grad_phi
std::function< const typename OutputTools< OutputShape >::VariablePhiDivergence &(const Assembly &, libMesh::FEType)> _div_phi_assembly_method
const FieldVariablePhiValue * _phi
unsigned int THREAD_ID
Definition: MooseTypes.h:209
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
uint8_t dof_id_type