https://mooseframework.inl.gov
MooseVariableData.h
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 #pragma once
11 
12 #include "MooseArray.h"
13 #include "MooseTypes.h"
14 #include "MooseVariableDataBase.h"
15 #include "Conversion.h"
16 
17 #include "libmesh/tensor_tools.h"
18 #include "libmesh/vector_value.h"
19 #include "libmesh/tensor_value.h"
20 #include "libmesh/type_n_tensor.h"
21 #include "libmesh/enum_fe_family.h"
22 #include "libmesh/fe_type.h"
23 #include "ADUtils.h"
24 
25 #include <functional>
26 #include <vector>
27 
28 namespace libMesh
29 {
30 class QBase;
31 class DofMap;
32 template <typename>
33 class DenseVector;
34 }
35 
37 using libMesh::DofMap;
38 using libMesh::QBase;
40 
41 class TimeIntegrator;
42 class Assembly;
43 class SubProblem;
44 template <typename>
45 class MooseVariableFE;
46 class SystemBase;
47 
48 template <typename OutputType>
49 class MooseVariableData : public MooseVariableDataBase<OutputType>
50 {
51 public:
52  // type for gradient, second and divergence of template class OutputType
56 
57  // shortcut for types storing values on quadrature points
63 
64  // shape function type for the template class OutputType
66 
67  // type for gradient, second and divergence of shape functions of template class OutputType
71 
72  // shortcut for types storing shape function values on quadrature points
78 
79  // shortcut for types storing test function values on quadrature points
80  // Note: here we assume the types are the same as of shape functions.
86 
87  // DoF value type for the template class OutputType
92 
94  SystemBase & sys,
95  THREAD_ID tid,
96  Moose::ElementType element_type,
97  const QBase * const & qrule_in,
98  const QBase * const & qrule_face_in,
99  const Node * const & node,
100  const Elem * const & elem);
101 
105  bool needsAD() const { return _need_ad; }
106 
111  void setGeometry(Moose::GeometryType gm_type);
112 
114 
118  void computeValues();
119 
124 
128  template <bool constant_monomial>
129  void computeAD(const unsigned int num_dofs, const unsigned int nqp);
130 
134  void computeNodalValues();
135 
137 
141  const FieldVariablePhiValue & phi() const { return *_phi; }
142 
146  const FieldVariablePhiValue & phiFace() const { return *_phi_face; }
147 
151  const FieldVariablePhiGradient & gradPhi() const { return *_grad_phi; }
152 
157  {
158  mooseAssert(var().fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY, "Not an array variable");
159  return _mapped_grad_phi;
160  }
161 
166 
171  {
172  mooseAssert(var().fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY, "Not an array variable");
173  return _mapped_grad_phi_face;
174  }
175 
179  const FieldVariablePhiSecond & secondPhi() const;
180 
184  const FieldVariablePhiSecond & secondPhiFace() const;
185 
189  const FieldVariablePhiCurl & curlPhi() const;
190 
194  const FieldVariablePhiCurl & curlPhiFace() const;
195 
199  const FieldVariablePhiDivergence & divPhi() const;
200 
204  const FieldVariablePhiDivergence & divPhiFace() const;
205 
210 
215 
219  std::size_t phiSize() const { return _phi->size(); }
220 
224  std::size_t phiFaceSize() const { return _phi_face->size(); }
225 
229  bool usesSecondPhi() const
230  {
232  }
233 
237  bool computingCurl() const { return _need_curl || _need_curl_old; }
238 
242  bool computingDiv() const { return _need_div || _need_div_old; }
243 
245 
246  bool isNodal() const override { return _is_nodal; }
247  bool hasDoFsOnNodes() const override { return _continuity != libMesh::DISCONTINUOUS; }
248  libMesh::FEContinuity getContinuity() const override { return _continuity; };
249  const Node * const & node() const { return _node; }
250  const dof_id_type & nodalDofIndex() const { return _nodal_dof_index; }
251  bool isNodalDefined() const { return _has_dof_indices; }
252 
256  const Elem * const & currentElem() const { return _elem; }
257 
261  const unsigned int & currentSide() const { return _current_side; }
262 
266  void prepareIC();
267 
269 
273  const FieldVariableGradient & gradSlnDot() const;
274 
278  const FieldVariableGradient & gradSlnDotDot() const;
279 
285 
290  const FieldVariableCurl & curlSln(Moose::SolutionState state) const;
291 
297 
299  {
300  _need_ad = _need_ad_u = true;
301  return _ad_u;
302  }
303 
305  {
306  _need_ad = _need_ad_grad_u = true;
307  return _ad_grad_u;
308  }
309 
311  {
312  _need_ad = _need_ad_grad_u_dot = true;
313 
314  if (!_time_integrator)
315  // If we don't have a time integrator (this will be the case for variables that are a part of
316  // the AuxiliarySystem) then we have no way to calculate _ad_grad_u_dot and we are just going
317  // to copy the values from _grad_u_dot. Of course in order to be able to do that we need to
318  // calculate _grad_u_dot
319  _need_grad_dot = true;
320 
321  return _ad_grad_u_dot;
322  }
323 
325  {
326  _need_ad = _need_ad_second_u = true;
327  secondPhi();
328  secondPhiFace();
329  return _ad_second_u;
330  }
331 
333 
335 
336  const FieldVariableValue & uDot() const;
337 
338  const FieldVariableValue & uDotDot() const;
339 
340  const FieldVariableValue & uDotOld() const;
341 
342  const FieldVariableValue & uDotDotOld() const;
343 
344  const VariableValue & duDotDu() const
345  {
346  _need_du_dot_du = true;
347  return _du_dot_du;
348  }
349 
350  const VariableValue & duDotDotDu() const
351  {
352  _need_du_dotdot_du = true;
353  return _du_dotdot_du;
354  }
355 
357  {
358  _need_ad = _need_ad_curl_u = true;
359  curlPhi();
360  curlPhiFace();
361  return _ad_curl_u;
362  }
363 
365 
366  const OutputType & nodalValueDot() const;
367  const OutputType & nodalValueDotDot() const;
368  const OutputType & nodalValueDotOld() const;
369  const OutputType & nodalValueDotDotOld() const;
370  const OutputType & nodalValueDuDotDu() const;
371  const OutputType & nodalValueDuDotDotDu() const;
372 
373  const typename Moose::ADType<OutputType>::type & adNodalValue() const;
374 
378  void setDofValues(const DenseVector<DofValue> & values);
379 
381 
384  void setDofValue(const DofValue & value, unsigned int index);
386 
391  DofValue getNodalValue(const Node & node, Moose::SolutionState state) const;
392  DofValue
393  getElementalValue(const Elem * elem, Moose::SolutionState state, unsigned int idx = 0) const;
394 
396 
397  void getDofIndices(const Elem * elem, std::vector<dof_id_type> & dof_indices) const;
398  const std::vector<dof_id_type> & dofIndices() const { return _dof_indices; }
399  unsigned int numberOfDofs() const { return _dof_indices.size(); }
400  void clearDofIndices() { _dof_indices.clear(); }
401 
405  void prepare();
406 
410  void reinitNode();
411 
415  void reinitAux();
416 
421  void reinitNodes(const std::vector<dof_id_type> & nodes);
422 
427  const DenseVector<libMesh::Number> & v) const;
428 
430 
431  const DofValues & dofValuesDot() const;
432  const DofValues & dofValuesDotOld() const;
433  const DofValues & dofValuesDotDot() const;
434  const DofValues & dofValuesDotDotOld() const;
437 
441  const ADDofValues & adDofValues() const;
442 
446  const ADDofValues & adDofValuesDot() const;
447 
449 
454  const FieldVariableValue & increment() const { return _increment; }
455 
460 
465 
466 private:
470  void fetchADDofValues();
471 
477  void assignADNodalValue();
478 
485  template <bool constant_monomial>
486  void computeValuesInternal();
487 
488  template <bool constant_monomial,
489  typename DestinationType,
490  typename ShapeType,
491  typename DofValuesType>
492  void fill(DestinationType & dest,
493  const ShapeType & phi,
494  const DofValuesType & dof_values,
495  unsigned int nqp,
496  std::size_t num_shapes);
497 
500 
502 
503  const unsigned int _var_num;
504 
506 
509 
511  bool _is_nodal;
512 
515 
518 
521 
524 
527 
529  mutable bool _need_second = false;
530  mutable bool _need_second_old = false;
531  mutable bool _need_second_older = false;
532  mutable bool _need_second_previous_nl = false;
533 
535  mutable bool _need_curl = false;
536  mutable bool _need_curl_old = false;
537  mutable bool _need_curl_older = false;
538 
540  mutable bool _need_div = false;
541  mutable bool _need_div_old = false;
542  mutable bool _need_div_older = false;
543 
545  mutable bool _need_ad = false;
546  mutable bool _need_ad_u = false;
547  mutable bool _need_ad_grad_u = false;
548  mutable bool _need_ad_grad_u_dot = false;
549  mutable bool _need_ad_second_u = false;
550  mutable bool _need_ad_curl_u = false;
551  mutable bool _need_ad_div_u = false;
552  mutable bool _need_ad_u_dot = false;
553  mutable bool _need_ad_u_dotdot = false;
554 
556 
560 
566 
571 
576 
590 
591  // time derivatives
592 
595 
598 
601 
604 
607 
610 
614  const QBase * const & _qrule;
615  const QBase * const & _qrule_face;
616 
617  // Shape function values, gradients, second derivatives
623 
624  // Mapped array phi
629 
630  // Values, gradients and second derivatives of shape function on faces
636 
639 
647 
648  // dual mortar
649  const bool _use_dual;
650 
657 
664 
671 
678 
682 
686 
687  std::function<const ADTemplateVariablePhiGradient<OutputShape> &(const Assembly &,
690  std::function<const ADTemplateVariablePhiGradient<OutputShape> &(const Assembly &,
693 
696 
700  const Node * const & _node;
701 
705  const Elem * const & _elem;
706 
708  const bool _displaced;
709 
711  const unsigned int & _current_side;
712 
715 
766 };
767 
769 
770 template <typename OutputType>
773 {
774  _need_ad = true;
775  return _ad_dof_values;
776 }
777 
778 template <typename OutputType>
781 {
782  _need_ad = _need_ad_u_dot = true;
783  if (!_time_integrator)
784  // See explanation in adUDot() body
785  _need_u_dot = true;
786  return _ad_dofs_dot;
787 }
788 
789 template <typename OutputType>
790 const typename Moose::ADType<OutputType>::type &
792 {
793  _need_ad = true;
794  return _ad_nodal_value;
795 }
796 
797 template <typename OutputType>
800 {
801  _need_ad = _need_ad_u_dot = true;
802 
803  if (!_time_integrator)
804  // If we don't have a time integrator (this will be the case for variables that are a part of
805  // the AuxiliarySystem) then we have no way to calculate _ad_u_dot and we are just going to
806  // copy the values from _u_dot. Of course in order to be able to do that we need to calculate
807  // _u_dot
808  _need_u_dot = true;
809 
810  return _ad_u_dot;
811 }
812 
813 template <typename OutputType>
816 {
817  _need_ad = _need_ad_u_dotdot = true;
818 
819  if (!_time_integrator)
820  // If we don't have a time integrator (this will be the case for variables that are a part
821  // of the AuxiliarySystem) then we have no way to calculate _ad_u_dotdot and we are just
822  // going to copy the values from _u_dotdot. Of course in order to be able to do that we need
823  // to calculate _u_dotdot
824  _need_u_dotdot = true;
825 
826  return _ad_u_dotdot;
827 }
828 
829 template <typename OutputType>
832 {
833  if (_element_type == Moose::ElementType::Neighbor || _element_type == Moose::ElementType::Lower)
834  mooseError("Unsupported element type: ", Moose::stringify(_element_type));
835  mooseAssert(_ad_grad_phi, "this should be non-null");
836  return *_ad_grad_phi;
837 }
838 
839 template <typename OutputType>
842 {
843  if (_element_type == Moose::ElementType::Neighbor || _element_type == Moose::ElementType::Lower)
844  mooseError("Unsupported element type: ", Moose::stringify(_element_type));
845  mooseAssert(_ad_grad_phi_face, "this should be non-null");
846  return *_ad_grad_phi_face;
847 }
const FieldVariableDivergence & divSln(Moose::SolutionState state) const
Local solution divergence getter.
typename OutputTools< typename Moose::ADType< T >::type >::VariableSecond ADTemplateVariableSecond
Definition: MooseTypes.h:654
FieldVariableCurl _curl_u_old
const Assembly & _assembly
bool _need_grad_dot
gradient dot flags
const ADTemplateVariableValue< OutputType > & adUDot() const
const FieldVariablePhiDivergence * _div_phi_face
typename OutputTools< typename Moose::ADType< T >::type >::VariableCurl ADTemplateVariableCurl
Definition: MooseTypes.h:656
const ADTemplateVariableSecond< OutputType > & adSecondSln() const
const FieldVariableValue & uDotDotOld() const
const OutputType & nodalValueDuDotDotDu() const
bool usesSecondPhi() const
Whether or not this variable is computing any second derivatives.
const MappedArrayVariablePhiGradient & arrayGradPhiFace() const
mapped_grad_phi_face getter
const ADDofValues & adDofValues() const
Return the AD dof values.
const FieldVariablePhiSecond * _second_phi
MooseArray< DofValue > DofValues
typename OutputTools< typename Moose::ADType< T >::type >::VariablePhiGradient ADTemplateVariablePhiGradient
Definition: MooseTypes.h:682
const unsigned int & currentSide() const
The current side.
MooseArray< std::vector< OutputShape > > FieldVariableTestCurl
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.
Keeps track of stuff related to assembling.
Definition: Assembly.h:109
Class for stuff related to variables.
Definition: Adaptivity.h:31
void prepare()
Get the dof indices corresponding to the current element.
MooseArray< OutputDivergence > FieldVariableDivergence
MappedArrayVariablePhiGradient _mapped_grad_phi
ADTemplateVariableValue< OutputType > _ad_u
AD u.
const MappedArrayVariablePhiGradient & arrayGradPhi() const
mapped_grad_phi getter
libMesh::TensorTools::DecrementRank< OutputShape >::type OutputShapeDivergence
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
MooseArray< std::vector< OutputShape > > FieldVariablePhiCurl
MooseArray< std::vector< OutputShapeDivergence > > FieldVariableTestDivergence
const ADTemplateVariableGradient< OutputType > & adGradSlnDot() const
OutputType type
Definition: MooseTypes.h:286
const Node *const & _node
The current node.
void computeConstantMonomialValues()
compute the values for const monomial variables
bool computingCurl() const
Whether or not this variable is computing the curl.
void computeValues()
compute the variable values
ADTemplateVariableSecond< OutputType > _ad_second_u
FieldVariableSecond _second_u_older
FieldVariableSecond _second_u_previous_nl
const FieldVariableValue & uDotDot() const
void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const
const DofValues & dofValuesDot() const
std::function< const typename OutputTools< OutputShape >::VariablePhiCurl &(const Assembly &, libMesh::FEType)> _curl_phi_assembly_method
bool _need_second
SolutionState second_u flags.
const OutputType & nodalValueDotOld() const
MooseArray< OutputGradient > FieldVariableGradient
const FieldVariablePhiGradient * _current_grad_phi
MooseArray< OutputSecond > FieldVariableSecond
const FieldVariableValue & increment() const
Increment getter.
MappedArrayVariablePhiGradient _mapped_grad_phi_face
std::function< const typename OutputTools< OutputShape >::VariablePhiCurl &(const Assembly &, libMesh::FEType)> _curl_phi_face_assembly_method
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:648
void computeIncrementAtNode(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change at the current node based on increment_vec.
const FieldVariablePhiDivergence * _div_phi
ADTemplateVariableGradient< OutputType > _ad_grad_u
dof_id_type _nodal_dof_index
The dof index for the current node.
const FieldVariablePhiCurl * _current_curl_phi
const FieldVariableValue & uDot() const
OutputTools< Real >::VariablePhiValue VariablePhiValue
Definition: MooseTypes.h:348
void assignADNodalValue()
Helper method for assigning nodal values from their corresponding solution values (dof values as they...
const ADTemplateVariableValue< OutputType > & adUDotDot() const
ADReal _ad_real_dummy
A dummy ADReal variable.
const FieldVariablePhiDivergence * _current_div_phi
const FieldVariablePhiCurl * _curl_phi_face
void computeNodalValues()
compute nodal things
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const
MooseArray< std::vector< OutputShapeSecond > > FieldVariablePhiSecond
const QBase * _current_qrule
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
ElementType
Definition: MooseTypes.h:812
const FieldVariablePhiCurl * _curl_phi
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
MooseArray< std::vector< OutputShapeSecond > > FieldVariableTestSecond
bool hasDoFsOnNodes() const override
Whether this data is associated with a variable that has DoFs on nodes.
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:42
ADTemplateVariableGradient< OutputType > _ad_grad_u_dot
FieldVariableValue _u_dot_old
u_dot_old (time derivative)
const unsigned int _var_num
void fill(DestinationType &dest, const ShapeType &phi, const DofValuesType &dof_values, unsigned int nqp, std::size_t num_shapes)
void prepareIC()
prepare the initial condition
Moose::ADType< DofValue >::type ADDofValue
MappedArrayVariablePhiGradient _mapped_grad_phi_neighbor
OutputTools< Real >::VariablePhiSecond VariablePhiSecond
Definition: MooseTypes.h:350
libMesh::TensorTools::DecrementRank< OutputType >::type OutputDivergence
MooseArray< std::vector< OutputShapeDivergence > > FieldVariablePhiDivergence
FieldVariableDivergence _div_u_old
const FieldVariableGradient & gradSlnDot() const
Local time derivative of solution gradient getter.
FieldVariableDivergence _div_u_older
libMesh::TensorTools::IncrementRank< OutputShapeGradient >::type OutputShapeSecond
void setDofValues(const DenseVector< DofValue > &values)
Set local DOF values and evaluate the values on quadrature points.
VariableValue _du_dot_du
derivative of u_dot wrt u
std::function< const typename OutputTools< OutputShape >::VariablePhiValue &(const Assembly &, libMesh::FEType)> _phi_face_assembly_method
void computeAD(const unsigned int num_dofs, const unsigned int nqp)
compute AD things
MappedArrayVariablePhiGradient _mapped_grad_phi_face_neighbor
std::function< const typename OutputTools< OutputShape >::VariablePhiSecond &(const Assembly &, libMesh::FEType)> _second_phi_face_assembly_method
MooseVariableData(const MooseVariableFE< 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)
const FieldVariablePhiSecond * _second_phi_face
const MooseArray< libMesh::Number > & dofValuesDuDotDu() const
FieldVariableValue _u_dot
u_dot (time derivative)
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
const FieldVariablePhiValue & phi() const
phi getter
ADTemplateVariableValue< OutputType > _ad_u_dot
libMesh::FEContinuity getContinuity() const override
Return the variable continuity.
FieldVariableSecond _second_u
second_u
const FieldVariablePhiDivergence & divPhi() const
divergence_phi getter
const DofValues & dofValuesDotDotOld() const
MooseArray< ADRealEigenVector > _ad_dofs_dot_eigen
void computeValuesInternal()
Internal method for computeValues() and computeConstantMonomialValues()
const OutputType & nodalValueDot() const
unsigned int numberOfDofs() const
const Elem *const & currentElem() const
The current element.
const FieldVariablePhiSecond & secondPhiFace() const
second_phi_face getter
ADDofValues _ad_dofs_dot
const QBase *const & _qrule_face
void computeIncrementAtQps(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change in solution at QPs based on increment_vec.
const VariableValue & duDotDotDu() const
const FieldVariablePhiGradient & gradPhiFace() const
grad_phi_face getter
const FieldVariableGradient & gradSlnDotDot() const
Local second time derivative of solution gradient getter.
ShapeType
Users of this template class must specify the type of shape functions that will be used in the Jacobi...
const ADTemplateVariableValue< OutputType > & adSln() const
bool _need_curl
curl flags
virtual const MooseVariableField< OutputType > & var() const
Moose::DOFType< OutputType >::type DofValue
std::vector< std::vector< Eigen::Map< RealDIMValue > > > MappedArrayVariablePhiGradient
Definition: MooseTypes.h:386
std::function< const typename OutputTools< OutputShape >::VariablePhiDivergence &(const Assembly &, libMesh::FEType)> _div_phi_face_assembly_method
FieldVariableGradient _grad_u_dotdot
const OutputType & nodalValueDuDotDu() const
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:651
OutputTools< Real >::VariablePhiCurl VariablePhiCurl
Definition: MooseTypes.h:351
bool _need_div
divergence flags
void reinitNode()
Prepare degrees of freedom for the current node.
const DofValues & dofValuesDotOld() const
bool isNodalDefined() const
const QBase *const & _qrule
The current qrule.
std::function< const ADTemplateVariablePhiGradient< OutputShape > &(const Assembly &, libMesh::FEType)> _ad_grad_phi_face_assembly_method
const TimeIntegrator * _time_integrator
Pointer to time integrator.
const MooseVariableFE< OutputType > & _var
A const reference to the owning MooseVariableFE object.
MooseArray< OutputType > FieldVariableCurl
FieldVariableValue _increment
Increment in the variable used in dampers.
FieldVariableDivergence _div_u
divergence_u
const Elem *const & _elem
The current elem.
std::string stringify(const T &t)
conversion to string
Definition: Conversion.h:64
const std::vector< dof_id_type > & dofIndices() const
const ADTemplateVariableCurl< OutputType > & adCurlSln() const
MooseArray< OutputType > FieldVariableValue
bool _need_ad
AD flags.
const ADTemplateVariablePhiGradient< OutputShape > & adGradPhi() const
ad_grad_phi getter
libMesh::FEContinuity _continuity
Continuity type of the variable.
OutputTools< Real >::VariablePhiDivergence VariablePhiDivergence
Definition: MooseTypes.h:352
std::size_t phiSize() const
Return phi size.
ADDofValues _ad_dofs_dotdot
const ADTemplateVariableGradient< OutputType > & adGradSln() const
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.
const Moose::ADType< OutputType >::type & adNodalValue() const
std::function< const typename OutputTools< OutputShape >::VariablePhiGradient &(const Assembly &, libMesh::FEType)> _grad_phi_assembly_method
const ADTemplateVariablePhiGradient< OutputShape > * _ad_grad_phi
Moose::ADType< OutputType >::type _ad_nodal_value
AD nodal value.
const dof_id_type & nodalDofIndex() const
const VariableValue & duDotDu() const
const OutputType & nodalValueDotDotOld() const
FieldVariableValue _u_dotdot
u_dotdot (second time derivative)
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:343
std::function< const typename OutputTools< OutputShape >::VariablePhiGradient &(const Assembly &, libMesh::FEType)> _grad_phi_face_assembly_method
const FieldVariablePhiValue * _phi_face
OutputTools< Real >::VariablePhiGradient VariablePhiGradient
Definition: MooseTypes.h:349
FieldVariableGradient _grad_u_dot
grad_u dots
const FieldVariablePhiValue & phiFace() const
phi_face getter
std::function< const typename OutputTools< OutputShape >::VariablePhiSecond &(const Assembly &, libMesh::FEType)> _second_phi_assembly_method
void fetchADDofValues()
Helper method for assigning the ad_dof* arrays.
const Node *const & node() const
FieldVariableValue _u_dotdot_old
u_dotdot_old (second time derivative)
Base class for time integrators.
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:78
ADDofValues _ad_dof_values
GeometryType
Definition: MooseTypes.h:277
DofValue getNodalValue(const Node &node, Moose::SolutionState state) const
const FieldVariableCurl & curlSln(Moose::SolutionState state) const
Local solution curl getter.
const DofValues & dofValuesDotDot() const
const FieldVariableValue & uDotOld() const
SolutionState
Definition: MooseTypes.h:261
const libMesh::FEType & _fe_type
FieldVariableSecond _second_u_old
Moose::ShapeType< OutputType >::type OutputShape
const ADTemplateVariablePhiGradient< OutputShape > * _ad_grad_phi_face
const FieldVariablePhiSecond & secondPhi() const
second_phi getter
const OutputType & nodalValueDotDot() const
bool computingDiv() const
Whether or not this variable is computing the divergence.
libMesh::TensorTools::IncrementRank< OutputShape >::type OutputShapeGradient
MooseArray< ADDofValue > ADDofValues
MooseArray< std::vector< OutputShape > > FieldVariablePhiValue
const bool _displaced
Whether this variable is being calculated on a displaced system.
const FieldVariablePhiSecond * _current_second_phi
const unsigned int & _current_side
The current element side.
void setDofValue(const DofValue &value, unsigned int index)
dof value setters
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.
const FieldVariablePhiCurl & curlPhiFace() const
curl_phi_face getter
const ADTemplateVariablePhiGradient< OutputShape > & adGradPhiFace() const
ad_grad_phi_face getter
const FieldVariablePhiGradient & gradPhi() const
grad_phi getter
MooseArray< std::vector< OutputShapeGradient > > FieldVariableTestGradient
VariableValue _du_dotdot_du
derivative of u_dotdot wrt u
MooseArray< ADRealEigenVector > _ad_dofs_dotdot_eigen
const FieldVariablePhiValue * _current_phi
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.
std::vector< dof_id_type > _dof_indices
The dof indices for the current element.
const ADTemplateVariablePhiGradient< OutputShape > * _current_ad_grad_phi
ADTemplateVariableValue< OutputType > _ad_u_dotdot
std::function< const typename OutputTools< OutputType >::VariablePhiValue &(const Assembly &, libMesh::FEType)> _phi_assembly_method
bool isNodal() const override
const FieldVariablePhiDivergence & divPhiFace() const
divergence_phi_face getter
const FieldVariableSecond & secondSln(Moose::SolutionState state) const
Local solution second spatial derivative getter.
FieldVariableCurl _curl_u_older
ADTemplateVariableCurl< OutputType > _ad_curl_u
std::size_t phiFaceSize() const
Return phiFace size.
bool _is_nodal
if variable is nodal
const FieldVariablePhiGradient * _grad_phi
void insertNodalValue(libMesh::NumericVector< libMesh::Number > &residual, const DofValue &v)
Write a nodal value to the passed-in solution vector.
std::function< const typename OutputTools< OutputShape >::VariablePhiDivergence &(const Assembly &, libMesh::FEType)> _div_phi_assembly_method
bool needsAD() const
Returns whether this data structure needs automatic differentiation calculations. ...
const FieldVariablePhiValue * _phi
const ADDofValues & adDofValuesDot() const
Return the AD time derivative values of degrees of freedom.
MooseArray< std::vector< OutputShape > > FieldVariableTestValue
unsigned int THREAD_ID
Definition: MooseTypes.h:237
MooseArray< std::vector< OutputShapeGradient > > FieldVariablePhiGradient
uint8_t dof_id_type
DofValue getElementalValue(const Elem *elem, Moose::SolutionState state, unsigned int idx=0) const
FieldVariableCurl _curl_u
curl_u