https://mooseframework.inl.gov
MooseVariableFE.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 "MooseTypes.h"
13 #include "MooseVariableFEBase.h"
14 #include "SubProblem.h"
15 #include "MooseMesh.h"
16 #include "MooseVariableField.h"
17 #include "MooseVariableData.h"
18 
19 #include "libmesh/numeric_vector.h"
20 #include "libmesh/dof_map.h"
21 #include "libmesh/elem.h"
22 #include "libmesh/quadrature.h"
23 #include "libmesh/dense_vector.h"
24 #include "libmesh/enum_fe_family.h"
25 
26 class TimeIntegrator;
27 template <typename>
32 
45 template <typename OutputType>
46 class MooseVariableFE : public MooseVariableField<OutputType>
47 {
48 public:
52 
58 
63 
71 
79 
84 
89 
91 
93 
94  void clearDofIndices() override;
95 
96  void prepare() override;
97  void prepareNeighbor() override;
98  void prepareLowerD() override;
99  virtual void prepareIC() override;
100 
101  void prepareAux() override;
102 
103  void reinitNode() override;
104  void reinitAux() override;
105  void reinitAuxNeighbor() override;
106 
107  void reinitNodes(const std::vector<dof_id_type> & nodes) override;
108  void reinitNodesNeighbor(const std::vector<dof_id_type> & nodes) override;
109 
115  bool usesPhi() const { return true; }
121  bool usesGradPhi() const { return true; }
122 
126  bool usesSecondPhi() const;
127 
132  bool usesSecondPhiNeighbor() const override final;
133 
137  bool computingSecond() const override final { return usesSecondPhi(); }
138 
142  bool computingCurl() const override final;
143 
147  bool computingDiv() const override final;
148 
149  bool isNodal() const override { return _element_data->isNodal(); }
150  bool hasDoFsOnNodes() const override { return _element_data->hasDoFsOnNodes(); }
151  libMesh::FEContinuity getContinuity() const override { return _element_data->getContinuity(); };
152  const Node * const & node() const { return _element_data->node(); }
153  const dof_id_type & nodalDofIndex() const override { return _element_data->nodalDofIndex(); }
154  virtual bool isNodalDefined() const override;
155 
156  const Node * const & nodeNeighbor() const { return _neighbor_data->node(); }
157  const dof_id_type & nodalDofIndexNeighbor() const override
158  {
159  return _neighbor_data->nodalDofIndex();
160  }
161  bool isNodalNeighborDefined() const;
162 
163  const Elem * const & currentElem() const override { return _element_data->currentElem(); }
164 
168  const unsigned int & currentSide() const { return _element_data->currentSide(); }
169 
173  const Elem * const & neighbor() const { return _neighbor_data->currentElem(); }
174 
175  virtual void getDofIndices(const Elem * elem,
176  std::vector<dof_id_type> & dof_indices) const override;
177  const std::vector<dof_id_type> & dofIndices() const final { return _element_data->dofIndices(); }
178  unsigned int numberOfDofs() const final { return _element_data->numberOfDofs(); }
179  const std::vector<dof_id_type> & dofIndicesNeighbor() const final
180  {
181  return _neighbor_data->dofIndices();
182  }
183  const std::vector<dof_id_type> & dofIndicesLower() const final
184  {
185  return _lower_data->dofIndices();
186  }
187 
188  void clearAllDofIndices() final;
189 
190  unsigned int numberOfDofsNeighbor() override { return _neighbor_data->dofIndices().size(); }
191 
192  const FieldVariablePhiValue & phi() const override { return _element_data->phi(); }
193  const FieldVariablePhiGradient & gradPhi() const override final
194  {
195  return _element_data->gradPhi();
196  }
198  {
199  return _element_data->arrayGradPhi();
200  }
201  const FieldVariablePhiSecond & secondPhi() const override final;
202  const FieldVariablePhiCurl & curlPhi() const override final;
203  const FieldVariablePhiDivergence & divPhi() const override final;
204 
205  const FieldVariablePhiValue & phiFace() const override final { return _element_data->phiFace(); }
206  const FieldVariablePhiGradient & gradPhiFace() const override final
207  {
208  return _element_data->gradPhiFace();
209  }
211  {
212  return _element_data->arrayGradPhiFace();
213  }
214  const FieldVariablePhiSecond & secondPhiFace() const override final;
215  const FieldVariablePhiCurl & curlPhiFace() const;
216  const FieldVariablePhiDivergence & divPhiFace() const;
217 
218  const FieldVariablePhiValue & phiNeighbor() const override final { return _neighbor_data->phi(); }
219  const FieldVariablePhiGradient & gradPhiNeighbor() const override final
220  {
221  return _neighbor_data->gradPhi();
222  }
224  {
225  return _neighbor_data->arrayGradPhi();
226  }
227  const FieldVariablePhiSecond & secondPhiNeighbor() const override final;
228  const FieldVariablePhiCurl & curlPhiNeighbor() const;
230 
231  const FieldVariablePhiValue & phiFaceNeighbor() const override final
232  {
233  return _neighbor_data->phiFace();
234  }
235  const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
236  {
237  return _neighbor_data->gradPhiFace();
238  }
240  {
241  return _neighbor_data->arrayGradPhiFace();
242  }
243  const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final;
246 
247  virtual const FieldVariablePhiValue & phiLower() const override { return _lower_data->phi(); }
248  const FieldVariablePhiGradient & gradPhiLower() const { return _lower_data->gradPhi(); }
249 
251  {
252  return _element_data->adGradPhi();
253  }
254 
256  {
257  return _element_data->adGradPhiFace();
258  }
259 
261  {
262  return _neighbor_data->adGradPhiFace();
263  }
264 
265  // damping
266  const FieldVariableValue & increment() const { return _element_data->increment(); }
267 
268  const FieldVariableValue & vectorTagValue(TagID tag) const override
269  {
270  return _element_data->vectorTagValue(tag);
271  }
273  {
274  return _element_data->vectorTagGradient(tag);
275  }
276  const DofValues & vectorTagDofValue(TagID tag) const override
277  {
278  return _element_data->vectorTagDofValue(tag);
279  }
280  const FieldVariableValue & matrixTagValue(TagID tag) const override
281  {
282  return _element_data->matrixTagValue(tag);
283  }
284 
286  const FieldVariableValue & sln() const override { return _element_data->sln(Moose::Current); }
287  const FieldVariableValue & slnOld() const override { return _element_data->sln(Moose::Old); }
288  const FieldVariableValue & slnOlder() const override { return _element_data->sln(Moose::Older); }
290 
292  const FieldVariableGradient & gradSln() const override
293  {
294  return _element_data->gradSln(Moose::Current);
295  }
296  const FieldVariableGradient & gradSlnOld() const override
297  {
298  return _element_data->gradSln(Moose::Old);
299  }
301  {
302  return _element_data->gradSln(Moose::Older);
303  }
305  {
306  return _element_data->gradSln(Moose::PreviousNL);
307  }
308 
310  const FieldVariableGradient & gradSlnDot() const { return _element_data->gradSlnDot(); }
311  const FieldVariableGradient & gradSlnDotDot() const { return _element_data->gradSlnDotDot(); }
312 
314  const FieldVariableSecond & secondSln() const { return _element_data->secondSln(Moose::Current); }
315  const FieldVariableSecond & secondSlnOld() const { return _element_data->secondSln(Moose::Old); }
317  {
318  return _element_data->secondSln(Moose::Older);
319  }
321  {
322  return _element_data->secondSln(Moose::PreviousNL);
323  }
324 
326  const FieldVariableCurl & curlSln() const { return _element_data->curlSln(Moose::Current); }
327  const FieldVariableCurl & curlSlnOld() const { return _element_data->curlSln(Moose::Old); }
328  const FieldVariableCurl & curlSlnOlder() const { return _element_data->curlSln(Moose::Older); }
329 
331  const FieldVariableDivergence & divSln() const { return _element_data->divSln(Moose::Current); }
332  const FieldVariableDivergence & divSlnOld() const { return _element_data->divSln(Moose::Old); }
334  {
335  return _element_data->divSln(Moose::Older);
336  }
337 
340  {
341  return _element_data->adSln();
342  }
343 
345  {
346  return _element_data->adGradSln();
347  }
349  {
350  return _element_data->adSecondSln();
351  }
353  {
354  return _element_data->adUDot();
355  }
357  {
358  return _element_data->adUDotDot();
359  }
361  {
362  return _element_data->adGradSlnDot();
363  }
365  {
366  return _element_data->adCurlSln();
367  }
368 
371  {
372  return _neighbor_data->adSln();
373  }
375  {
376  return _neighbor_data->adGradSln();
377  }
379  {
380  return _neighbor_data->adSecondSln();
381  }
383  {
384  return _neighbor_data->adUDot();
385  }
387  {
388  return _neighbor_data->adUDotDot();
389  }
391  {
392  return _neighbor_data->adGradSlnDot();
393  }
395  {
396  return _neighbor_data->adCurlSln();
397  }
398 
400  const FieldVariableValue & uDot() const { return _element_data->uDot(); }
401  const FieldVariableValue & uDotDot() const { return _element_data->uDotDot(); }
402  const FieldVariableValue & uDotOld() const { return _element_data->uDotOld(); }
403  const FieldVariableValue & uDotDotOld() const { return _element_data->uDotDotOld(); }
404  const VariableValue & duDotDu() const { return _element_data->duDotDu(); }
405  const VariableValue & duDotDotDu() const { return _element_data->duDotDotDu(); }
406 
408  const FieldVariableValue & slnNeighbor() const override
409  {
410  return _neighbor_data->sln(Moose::Current);
411  }
412  const FieldVariableValue & slnOldNeighbor() const override
413  {
414  return _neighbor_data->sln(Moose::Old);
415  }
418  {
419  return _neighbor_data->sln(Moose::PreviousNL);
420  }
421 
423  const FieldVariableGradient & gradSlnNeighbor() const override
424  {
425  return _neighbor_data->gradSln(Moose::Current);
426  }
427  const FieldVariableGradient & gradSlnOldNeighbor() const override
428  {
429  return _neighbor_data->gradSln(Moose::Old);
430  }
432  {
433  return _neighbor_data->gradSln(Moose::Older);
434  }
436  {
437  return _neighbor_data->gradSln(Moose::PreviousNL);
438  }
439 
441  const FieldVariableGradient & gradSlnNeighborDot() const { return _neighbor_data->gradSlnDot(); }
443  {
444  return _neighbor_data->gradSlnDotDot();
445  }
446 
449  {
450  return _neighbor_data->secondSln(Moose::Current);
451  }
453  {
454  return _neighbor_data->secondSln(Moose::Old);
455  }
457  {
458  return _neighbor_data->secondSln(Moose::Older);
459  }
461  {
462  return _neighbor_data->secondSln(Moose::PreviousNL);
463  }
464 
467  {
468  return _neighbor_data->curlSln(Moose::Current);
469  }
471  {
472  return _neighbor_data->curlSln(Moose::Old);
473  }
475  {
476  return _neighbor_data->curlSln(Moose::Older);
477  }
478 
481  {
482  return _neighbor_data->divSln(Moose::Current);
483  }
485  {
486  return _neighbor_data->divSln(Moose::Old);
487  }
489  {
490  return _neighbor_data->divSln(Moose::Older);
491  }
492 
494  const FieldVariableValue & uDotNeighbor() const { return _neighbor_data->uDot(); }
495  const FieldVariableValue & uDotDotNeighbor() const { return _neighbor_data->uDotDot(); }
496  const FieldVariableValue & uDotOldNeighbor() const { return _neighbor_data->uDotOld(); }
497  const FieldVariableValue & uDotDotOldNeighbor() const { return _neighbor_data->uDotDotOld(); }
498  const VariableValue & duDotDuNeighbor() const { return _neighbor_data->duDotDu(); }
499  const VariableValue & duDotDotDuNeighbor() const { return _neighbor_data->duDotDotDu(); }
500 
502  const ADTemplateVariableValue<OutputType> & adSlnLower() const { return _lower_data->adSln(); }
503  const FieldVariableValue & slnLower() const { return _lower_data->sln(Moose::Current); }
504  const FieldVariableValue & slnLowerOld() const { return _lower_data->sln(Moose::Old); }
505 
507  virtual void computeElemValues() override;
508  virtual void computeElemValuesFace() override;
509  virtual void computeNeighborValuesFace() override;
510  virtual void computeNeighborValues() override;
511  virtual void computeLowerDValues() override;
512 
513  virtual void setNodalValue(const OutputType & value, unsigned int idx = 0) override;
514 
515  virtual void setDofValue(const DofValue & value, unsigned int index) override;
516 
520  virtual void setDofValues(const DenseVector<DofValue> & values) override;
521  virtual void setLowerDofValues(const DenseVector<DofValue> & values) override;
522 
527 
531  DofValue getNodalValue(const Node & node) const;
535  DofValue getNodalValueOld(const Node & node) const;
539  DofValue getNodalValueOlder(const Node & node) const;
546  DofValue getElementalValue(const Elem * elem, unsigned int idx = 0) const;
553  DofValue getElementalValueOld(const Elem * elem, unsigned int idx = 0) const;
560  DofValue getElementalValueOlder(const Elem * elem, unsigned int idx = 0) const;
561 
565  virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) override;
566  virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) override;
567 
571  virtual void add(libMesh::NumericVector<libMesh::Number> & vector) override;
572 
576  void addSolution(const DenseVector<libMesh::Number> & v);
577 
581  void addSolutionNeighbor(const DenseVector<libMesh::Number> & v);
582 
583  const DofValues & dofValue() const;
584  const DofValues & dofValues() const override;
585  const DofValues & dofValuesOld() const override;
586  const DofValues & dofValuesOlder() const override;
587  const DofValues & dofValuesPreviousNL() const override;
588  const DofValues & dofValuesNeighbor() const override;
589  const DofValues & dofValuesOldNeighbor() const override;
590  const DofValues & dofValuesOlderNeighbor() const override;
591  const DofValues & dofValuesPreviousNLNeighbor() const override;
592  const DofValues & dofValuesDot() const override;
593  const DofValues & dofValuesDotNeighbor() const override;
594  const DofValues & dofValuesDotNeighborResidual() const;
595  const DofValues & dofValuesDotOld() const override;
596  const DofValues & dofValuesDotOldNeighbor() const override;
597  const DofValues & dofValuesDotDot() const override;
598  const DofValues & dofValuesDotDotNeighbor() const override;
600  const DofValues & dofValuesDotDotOld() const override;
601  const DofValues & dofValuesDotDotOldNeighbor() const override;
602  const MooseArray<libMesh::Number> & dofValuesDuDotDu() const override;
603  const MooseArray<libMesh::Number> & dofValuesDuDotDuNeighbor() const override;
604  const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const override;
606 
607  const ADDofValues & adDofValues() const override;
608  const ADDofValues & adDofValuesNeighbor() const override;
609  const ADDofValues & adDofValuesDot() const override;
610 
615 
620 
627  OutputType getValue(const Elem * elem, const std::vector<std::vector<OutputShape>> & phi) const;
628 
636  const Elem * elem,
637  const std::vector<std::vector<typename OutputTools<OutputType>::OutputShapeGradient>> &
638  grad_phi) const;
639 
643  virtual std::size_t phiSize() const final { return _element_data->phiSize(); }
647  virtual std::size_t phiFaceSize() const final { return _element_data->phiFaceSize(); }
651  virtual std::size_t phiNeighborSize() const final { return _neighbor_data->phiSize(); }
655  virtual std::size_t phiFaceNeighborSize() const final { return _neighbor_data->phiFaceSize(); }
656 
657  std::size_t phiLowerSize() const final { return _lower_data->phiSize(); }
658 
662  const OutputType & nodalValue() const;
663  const OutputType & nodalValueOld() const;
664  const OutputType & nodalValueOlder() const;
665  const OutputType & nodalValuePreviousNL() const;
666  const OutputType & nodalValueDot() const;
667  const OutputType & nodalValueDotDot() const;
668  const OutputType & nodalValueDotOld() const;
669  const OutputType & nodalValueDotDotOld() const;
670  const OutputType & nodalValueDuDotDu() const;
671  const OutputType & nodalValueDuDotDotDu() const;
672  const OutputType & nodalValueNeighbor() const;
673  const OutputType & nodalValueOldNeighbor() const;
674  const OutputType & nodalValueOlderNeighbor() const;
675  const OutputType & nodalValuePreviousNLNeighbor() const;
676  const OutputType & nodalValueDotNeighbor() const;
677  const OutputType & nodalValueDotNeighborResidual() const;
678  const OutputType & nodalValueDotDotNeighbor() const;
679  const OutputType & nodalValueDotDotNeighborResidual() const;
680  const OutputType & nodalValueDotOldNeighbor() const;
681  const OutputType & nodalValueDotDotOldNeighbor() const;
682  const OutputType & nodalValueDuDotDuNeighbor() const;
683  const OutputType & nodalValueDuDotDotDuNeighbor() const;
684 
685  const MooseArray<OutputType> & nodalValueArray() const override
686  {
687  return _element_data->nodalValueArray(Moose::Current);
688  }
690  {
691  return _element_data->nodalValueArray(Moose::Old);
692  }
694  {
695  return _element_data->nodalValueArray(Moose::Older);
696  }
697 
698  const DofValues & nodalVectorTagValue(TagID tag) const override;
699  const DofValues & nodalMatrixTagValue(TagID tag) const override;
700 
701  const typename Moose::ADType<OutputType>::type & adNodalValue() const;
702 
703  virtual void computeNodalValues() override;
704  virtual void computeNodalNeighborValues() override;
705 
706  unsigned int oldestSolutionStateRequested() const override final;
707 
708  void setActiveTags(const std::set<TagID> & vtags) override;
709 
710  virtual void meshChanged() override;
711  virtual void residualSetup() override;
712  virtual void jacobianSetup() override;
713 
714  virtual void sizeMatrixTagData() override;
715 
716  bool supportsFaceArg() const override final { return true; }
717  bool supportsElemSideQpArg() const override final { return true; }
718 
719 protected:
721 
723  std::unique_ptr<MooseVariableData<OutputType>> _element_data;
724 
726  std::unique_ptr<MooseVariableData<OutputType>> _neighbor_data;
727 
729  std::unique_ptr<MooseVariableData<OutputType>> _lower_data;
730 
742 
749  ValueType
750  faceEvaluate(const FaceArg &, const StateArg &, const std::vector<ValueType> & cache_data) const;
751 
752  ValueType evaluate(const ElemQpArg & elem_qp, const StateArg & state) const override final;
753  ValueType evaluate(const ElemSideQpArg & elem_side_qp,
754  const StateArg & state) const override final;
755  ValueType evaluate(const ElemArg &, const StateArg &) const override final;
756  ValueType evaluate(const ElemPointArg &, const StateArg &) const override final;
757  ValueType evaluate(const NodeArg & node_arg, const StateArg & state) const override final;
758  ValueType evaluate(const FaceArg &, const StateArg &) const override final;
759 
760  GradientType evaluateGradient(const ElemQpArg & elem_qp, const StateArg & state) const override;
761  GradientType evaluateGradient(const ElemSideQpArg & elem_side_qp,
762  const StateArg & state) const override final;
763  GradientType evaluateGradient(const ElemArg &, const StateArg &) const override final;
764 
765  DotType evaluateDot(const ElemQpArg & elem_qp, const StateArg & state) const override final;
766  DotType evaluateDot(const ElemSideQpArg & elem_side_qp,
767  const StateArg & state) const override final;
768  DotType evaluateDot(const ElemArg &, const StateArg &) const override final;
769  DotType evaluateDot(const FaceArg &, const StateArg &) const override final;
770 
771  GradientType evaluateGradDot(const ElemArg &, const StateArg &) const override final;
772 
773 private:
778  template <typename Shapes, typename Solution, typename GradShapes, typename GradSolution>
779  void computeSolution(const Elem * elem,
780  unsigned int n_qp,
781  const StateArg & state,
782  const Shapes & phi,
783  Solution & local_soln,
784  const GradShapes & grad_phi,
785  GradSolution & grad_local_soln,
786  Solution & dot_local_soln,
787  GradSolution & grad_dot_local_soln) const;
788 
792  void
793  evaluateOnElement(const ElemQpArg & elem_qp, const StateArg & state, bool cache_eligible) const;
794 
798  void evaluateOnElementSide(const ElemSideQpArg & elem_side_qp,
799  const StateArg & state,
800  bool cache_eligible) const;
801 
805  mutable const Elem * _current_elem_qp_functor_elem = nullptr;
806 
808  mutable std::vector<ValueType> _current_elem_qp_functor_sln;
809 
811  mutable std::vector<GradientType> _current_elem_qp_functor_gradient;
812 
814  mutable std::vector<DotType> _current_elem_qp_functor_dot;
815 
817  mutable std::vector<GradientType> _current_elem_qp_functor_grad_dot;
818 
822  mutable std::pair<const Elem *, unsigned int> _current_elem_side_qp_functor_elem_side{
823  nullptr, libMesh::invalid_uint};
824 
826  mutable std::vector<ValueType> _current_elem_side_qp_functor_sln;
827 
829  mutable std::vector<GradientType> _current_elem_side_qp_functor_gradient;
830 
832  mutable std::vector<DotType> _current_elem_side_qp_functor_dot;
833 
836  mutable std::vector<GradientType> _current_elem_side_qp_functor_grad_dot;
837 };
838 
839 template <typename OutputType>
840 inline const typename MooseVariableFE<OutputType>::ADDofValues &
842 {
843  return _element_data->adDofValues();
844 }
845 
846 template <typename OutputType>
847 inline const typename MooseVariableFE<OutputType>::ADDofValues &
849 {
850  return _neighbor_data->adDofValues();
851 }
852 
853 template <typename OutputType>
854 inline const typename MooseVariableFE<OutputType>::ADDofValues &
856 {
857  return _element_data->adDofValuesDot();
858 }
859 
860 template <typename OutputType>
861 inline const typename Moose::ADType<OutputType>::type &
863 {
864  return _element_data->adNodalValue();
865 }
866 
867 template <typename OutputType>
868 void
869 MooseVariableFE<OutputType>::setActiveTags(const std::set<TagID> & vtags)
870 {
871  _element_data->setActiveTags(vtags);
872  _neighbor_data->setActiveTags(vtags);
873  _lower_data->setActiveTags(vtags);
874 }
875 
876 // Declare all the specializations, as the template specialization declarations below must know
877 template <>
879 template <>
881 template <>
883 template <>
886  const std::vector<std::vector<Real>> & phi) const;
887 template <>
889  const Elem * elem, const std::vector<std::vector<RealVectorValue>> & grad_phi) const;
890 template <>
892  const StateArg &,
893  bool) const;
894 template <>
896  const StateArg &,
897  bool) const;
898 template <>
900 MooseVariableFE<RealEigenVector>::evaluate(const ElemQpArg &, const StateArg &) const;
901 template <>
903 MooseVariableFE<RealEigenVector>::evaluate(const ElemSideQpArg &, const StateArg &) const;
904 template <>
906 MooseVariableFE<RealEigenVector>::evaluateGradient(const ElemQpArg &, const StateArg &) const;
907 template <>
909 MooseVariableFE<RealEigenVector>::evaluateGradient(const ElemSideQpArg &, const StateArg &) const;
910 template <>
912 MooseVariableFE<RealEigenVector>::evaluateDot(const ElemQpArg &, const StateArg &) const;
913 template <>
915 MooseVariableFE<RealEigenVector>::evaluateDot(const ElemSideQpArg &, const StateArg &) const;
916 
917 // Prevent implicit instantiation in other translation units where these classes are used
918 extern template class MooseVariableFE<Real>;
919 extern template class MooseVariableFE<RealVectorValue>;
920 extern template class MooseVariableFE<RealEigenVector>;
typename MooseVariableField< Real >::FieldVariablePhiValue FieldVariablePhiValue
typename OutputTools< typename Moose::ADType< T >::type >::VariableSecond ADTemplateVariableSecond
Definition: MooseTypes.h:654
const ADDofValues & adDofValuesDot() const override
Return the AD time derivatives at dofs.
const FieldVariableValue & vectorTagValue(TagID tag) const override
tag values getters
void computeSolution(const Elem *elem, unsigned int n_qp, const StateArg &state, const Shapes &phi, Solution &local_soln, const GradShapes &grad_phi, GradSolution &grad_local_soln, Solution &dot_local_soln, GradSolution &grad_dot_local_soln) const
Compute the solution, gradient, time derivative, and gradient of the time derivative with provided sh...
virtual void setDofValues(const DenseVector< DofValue > &values) override
Set local DOF values and evaluate the values on quadrature points.
const ADTemplateVariableGradient< OutputType > & adGradSlnDot() const override
AD grad of time derivative solution getter.
const FieldVariableValue & uDotDotNeighbor() const
const FieldVariableCurl & curlSlnNeighbor() const
neighbor solution curls
const Node *const & node() const
typename MooseVariableField< Real >::OutputShapeGradient OutputShapeGradient
typename OutputTools< typename Moose::ADType< T >::type >::VariableCurl ADTemplateVariableCurl
Definition: MooseTypes.h:656
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const override
const FieldVariableSecond & secondSlnOlderNeighbor() const
GradientType evaluateGradDot(const ElemArg &, const StateArg &) const override final
Evaluate the functor gradient-dot with a given element.
bool computingCurl() const override final
Whether or not this variable is computing the curl.
bool isNodalNeighborDefined() const
const OutputType & nodalValueDuDotDotDu() const
const DofValues & nodalVectorTagValue(TagID tag) const override
bool computingDiv() const override final
Whether or not this variable is computing the divergence.
const ADTemplateVariableTestGradient< OutputShape > & adGradPhiFace() const
OutputTools< OutputType >::OutputGradient getGradient(const Elem *elem, const std::vector< std::vector< typename OutputTools< OutputType >::OutputShapeGradient >> &grad_phi) const
Compute the variable gradient value at a point on an element.
ValueType evaluate(const ElemQpArg &elem_qp, const StateArg &state) const override final
const FieldVariableGradient & gradSlnPreviousNL() const
typename MooseVariableField< Real >::OutputGradient OutputGradient
const unsigned int invalid_uint
const FieldVariablePhiValue & phiNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element.
const OutputType & nodalValueOld() const
Class for stuff related to variables.
Definition: Adaptivity.h:31
typename MooseVariableField< Real >::FieldVariablePhiSecond FieldVariablePhiSecond
unsigned int TagID
Definition: MooseTypes.h:238
const OutputType & nodalValueDuDotDu() const
const OutputType & nodalValueDotDotOld() const
const FieldVariablePhiDivergence & divPhiNeighbor() const
const OutputType & nodalValueDot() const
Class for stuff related to variables.
const MooseArray< libMesh::Number > & dofValuesDuDotDu() const override
const FieldVariableSecond & secondSlnPreviousNLNeighbor() const
OutputType getValue(const Elem *elem, const std::vector< std::vector< OutputShape >> &phi) const
Compute the variable value at a point on an element.
const FieldVariableCurl & curlSlnOldNeighbor() const
DofValue getElementalValueOlder(const Elem *elem, unsigned int idx=0) const
Get the older value of this variable on an element.
virtual void setNodalValue(const OutputType &value, unsigned int idx=0) override
virtual void computeLowerDValues() override
compute values at quadrature points on the lower dimensional element
const FieldVariableGradient & gradSlnNeighborDotDot() const
std::vector< DotType > _current_elem_qp_functor_dot
The values of the time derivative for the _current_elem_qp_functor_elem.
const DofValues & dofValuesDotDotNeighborResidual() const
const OutputType & nodalValueNeighbor() const
typename MooseVariableField< Real >::FieldVariablePhiCurl FieldVariablePhiCurl
const FieldVariableGradient & gradSln() const override
element gradients
const FieldVariablePhiCurl & curlPhiFaceNeighbor() const
const FieldVariableValue & slnNeighbor() const override
neighbor solutions
const MooseArray< libMesh::Number > & dofValuesDuDotDuNeighbor() const override
const VariableValue & duDotDuNeighbor() const
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
void reinitAux() override
const FieldVariableSecond & secondSlnOld() const
const FieldVariableValue & slnOlderNeighbor() const
const OutputType & nodalValuePreviousNLNeighbor() const
void addSolution(const DenseVector< libMesh::Number > &v)
Add passed in local DOF values onto the current solution.
const ADTemplateVariableSecond< OutputType > & adSecondSlnNeighbor() const override
AD second neighbor solution getter.
virtual std::size_t phiFaceNeighborSize() const final
Return phiFaceNeighbor size.
void reinitNodes(const std::vector< dof_id_type > &nodes) override
const std::vector< dof_id_type > & dofIndicesLower() const final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
virtual void sizeMatrixTagData() override
Size data structures related to matrix tagging.
const OutputType & nodalValueDuDotDuNeighbor() const
const OutputType & nodalValueDuDotDotDuNeighbor() const
std::unique_ptr< MooseVariableData< OutputType > > _neighbor_data
Holder for all the data associated with the neighbor element.
virtual void setDofValue(const DofValue &value, unsigned int index) override
Degree of freedom value setters.
void prepareNeighbor() override
Prepare the neighbor element degrees of freedom.
MooseVariableFE(const InputParameters &parameters)
const FieldVariablePhiGradient & gradPhiNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element.
const OutputType & nodalValueDotDotNeighbor() const
const DofValues & dofValuesOlder() const override
void insertNodalValue(libMesh::NumericVector< libMesh::Number > &residual, const DofValue &v)
Write a nodal value to the passed-in solution vector.
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:648
const FieldVariableValue & slnLowerOld() const
const FieldVariableDivergence & divSlnOld() const
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const DofValues & dofValuesDotDotOldNeighbor() const override
const MooseArray< libMesh::Number > & dofValuesDuDotDotDuNeighbor() const override
const FieldVariableGradient & gradSlnPreviousNLNeighbor() const
typename MooseVariableField< Real >::FieldVariableValue FieldVariableValue
const OutputType & nodalValueDotDotOldNeighbor() const
const ADTemplateVariableValue< OutputType > & adSlnNeighbor() const override
neighbor AD
const FieldVariablePhiSecond & secondPhi() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s elemental shape functions.
const MappedArrayVariablePhiGradient & arrayGradPhiFaceNeighbor() const
const FieldVariableValue & uDotOld() const
typename MooseVariableField< Real >::OutputShapeSecond OutputShapeSecond
const FieldVariablePhiDivergence & divPhi() const override final
Divergence of the shape functions.
typename MooseVariableField< Real >::FieldVariableGradient FieldVariableGradient
bool hasDoFsOnNodes() const override
Does this variable have DoFs on nodes.
A structure that is used to evaluate Moose functors at an arbitrary physical point contained within a...
virtual void prepareIC() override
Prepare the initial condition.
bool usesSecondPhiNeighbor() const override final
Whether or not this variable is actually using the shape function second derivative on a neighbor...
const ADTemplateVariableValue< OutputType > & adUDotNeighbor() const override
AD neighbor time derivative getter.
const FieldVariableCurl & curlSlnOld() const
const ADTemplateVariableValue< OutputType > & adUDot() const override
AD time derivative getter.
const ADTemplateVariableGradient< OutputType > & adGradSlnNeighborDot() const override
AD grad of time derivative neighbor solution getter.
const DofValues & dofValuesOld() const override
typename MooseVariableField< Real >::OutputSecond OutputSecond
void addSolutionNeighbor(const DenseVector< libMesh::Number > &v)
Add passed in local neighbor DOF values onto the current solution.
const DofValues & nodalMatrixTagValue(TagID tag) const override
virtual void jacobianSetup() override
typename FunctorReturnType< T, FunctorEvaluationKind::Gradient >::type GradientType
This rigmarole makes it so that a user can create functors that return containers (std::vector...
Definition: MooseFunctor.h:149
std::vector< GradientType > _current_elem_qp_functor_grad_dot
The values of the gradient of the time derivative for the _current_elem_qp_functor_elem.
const ADTemplateVariableSecond< OutputType > & adSecondSln() const override
AD second solution getter.
bool supportsFaceArg() const override final
Whether this functor supports evaluation with FaceArg.
std::pair< const Elem *, unsigned int > _current_elem_side_qp_functor_elem_side
Keep track of the current elem-side-qp functor element and side in order to enable local caching (e...
const DofValues & dofValuesDotOldNeighbor() const override
const FieldVariableGradient & gradSlnOlderNeighbor() const
const ADDofValues & adDofValuesNeighbor() const override
Return the AD neighbor dof values.
virtual void computeNodalValues() override
Compute nodal values of this variable.
const DofValues & dofValuesDot() const override
const FieldVariableSecond & secondSlnOldNeighbor() const
bool usesSecondPhi() const
Whether or not this variable is computing any second derivatives.
const MooseArray< OutputType > & nodalValueOldArray() const override
const FieldVariablePhiCurl & curlPhiNeighbor() const
virtual void insertLower(libMesh::NumericVector< libMesh::Number > &vector) override
Insert the currently cached degree of freedom values for a lower-dimensional element into the provide...
const FieldVariableValue & uDot() const
element dots
const DofValues & dofValue() const
MooseVariableFE< RealVectorValue > VectorMooseVariable
const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element face.
GradientType evaluateGradient(const ElemQpArg &elem_qp, const StateArg &state) const override
void prepareLowerD() override
Prepare a lower dimensional element&#39;s degrees of freedom.
const FieldVariableCurl & curlSln() const
element curls
typename MooseVariableField< Real >::FieldVariableTestValue FieldVariableTestValue
std::vector< DotType > _current_elem_side_qp_functor_dot
The values of the time derivative for the _current_elem_side_qp_functor_elem_side.
ValueType faceEvaluate(const FaceArg &, const StateArg &, const std::vector< ValueType > &cache_data) const
A common method that both evaluate(FaceArg) and evaluateDot(FaceArg) can call.
unsigned int numberOfDofsNeighbor() override
const VariableValue & duDotDotDu() const
const ADTemplateVariableCurl< OutputType > & adCurlSln() const override
AD curl solution getter.
libMesh::TensorTools::DecrementRank< OutputShape >::type OutputShapeDivergence
const MappedArrayVariablePhiGradient & arrayGradPhi() const
typename MooseVariableField< Real >::FieldVariableDivergence FieldVariableDivergence
const FieldVariablePhiCurl & curlPhiFace() const
const MooseArray< OutputType > & nodalValueOlderArray() const override
const OutputType & nodalValueDotDot() const
const FieldVariableValue & uDotDotOldNeighbor() const
const OutputType & nodalValueDotDotNeighborResidual() const
DofValue getNodalValueOlder(const Node &node) const
Get the t-2 value of this variable at given node.
const DofValues & dofValuesOlderNeighbor() const override
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
const ADTemplateVariableTestGradient< OutputShape > & adGradPhiFaceNeighbor() const
const FieldVariableDivergence & divSlnOlder() const
const FieldVariableSecond & secondSlnPreviousNL() const
const OutputType & nodalValueOlder() const
A structure defining a "face" evaluation calling argument for Moose functors.
void reinitAuxNeighbor() override
DofValue getNodalValueOld(const Node &node) const
Get the old value of this variable at given node.
const MappedArrayVariablePhiGradient & arrayGradPhiNeighbor() const
const DofValues & dofValuesDotDotOld() const override
virtual void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const override
virtual void computeElemValuesFace() override
Compute values at facial quadrature points.
const FieldVariableValue & slnPreviousNL() const
const OutputType & nodalValueDotNeighbor() const
unsigned int numberOfDofs() const final
Get the number of local DoFs.
const FieldVariableValue & uDotDot() const
const FieldVariableValue & sln() const override
element solutions
void setActiveTags(const std::set< TagID > &vtags) override
Set the active vector tags.
const FieldVariableSecond & secondSlnNeighbor() const
neighbor solution seconds
const FieldVariablePhiValue & phi() const override
Return the variable&#39;s elemental shape functions.
libMesh::FEContinuity getContinuity() const override
Return the continuity of this variable.
const DofValues & dofValuesDotNeighborResidual() const
MooseVariableFE< RealEigenVector > ArrayMooseVariable
std::vector< std::vector< Eigen::Map< RealDIMValue > > > MappedArrayVariablePhiGradient
Definition: MooseTypes.h:386
std::size_t phiLowerSize() const final
Return the number of shape functions on the lower dimensional element for this variable.
typename OutputTools< typename Moose::ADType< T >::type >::VariableTestGradient ADTemplateVariableTestGradient
Definition: MooseTypes.h:676
typename MooseVariableField< Real >::FieldVariableTestDivergence FieldVariableTestDivergence
const FieldVariableSecond & secondSlnOlder() const
std::vector< GradientType > _current_elem_side_qp_functor_grad_dot
The values of the gradient of the time derivative for the _current_elem_side_qp_functor_elem_side.
typename MooseVariableField< Real >::FieldVariableTestSecond FieldVariableTestSecond
typename MooseVariableField< Real >::OutputShapeDivergence OutputShapeDivergence
const FieldVariableValue & uDotDotOld() const
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:651
const unsigned int & currentSide() const
Current side this variable is being evaluated on.
virtual void computeNodalNeighborValues() override
Compute nodal values of this variable in the neighbor.
typename MooseVariableField< Real >::FieldVariableSecond FieldVariableSecond
const FieldVariableGradient & vectorTagGradient(TagID tag) const
A structure that is used to evaluate Moose functors logically at an element/cell center.
virtual void meshChanged() override
Called on this object when the mesh changes.
const ADTemplateVariableValue< OutputType > & adSln() const override
AD.
std::vector< ValueType > _current_elem_qp_functor_sln
The values of the solution for the _current_elem_qp_functor_elem.
const FieldVariablePhiSecond & secondPhiFace() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s shape functions on an element face...
const Node *const & nodeNeighbor() const
Argument for requesting functor evaluation at a quadrature point location in an element.
const FieldVariableValue & uDotNeighbor() const
neighbor dots
const OutputType & nodalValue() const
Methods for retrieving values of variables at the nodes.
const FieldVariableValue & increment() const
typename MooseVariableField< Real >::FieldVariableCurl FieldVariableCurl
virtual void add(libMesh::NumericVector< libMesh::Number > &vector) override
Add the current local DOF values to the input vector.
void reinitNodesNeighbor(const std::vector< dof_id_type > &nodes) override
const DofValues & dofValuesNeighbor() const override
const FieldVariablePhiSecond & secondPhiNeighbor() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s shape functions on a neighboring ele...
std::unique_ptr< MooseVariableData< OutputType > > _lower_data
Holder for all the data associated with the lower dimeensional element.
const FieldVariableSecond & secondSln() const
element seconds
virtual std::size_t phiFaceSize() const final
Return phiFace size.
const DofValues & dofValuesDotDot() const override
typename MooseVariableField< Real >::FieldVariableTestGradient FieldVariableTestGradient
void computeIncrementAtQps(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change in solution at QPs based on increment_vec.
void clearDofIndices() override
Clear out the dof indices.
const FieldVariableValue & slnLower() const
const FieldVariableGradient & gradSlnOldNeighbor() const override
MooseVariableFE< Real > MooseVariable
const ADTemplateVariableValue< OutputType > & adUDotDot() const override
AD second time derivative getter.
const FieldVariableGradient & gradSlnOld() const override
typename MooseVariableField< Real >::FieldVariablePhiGradient FieldVariablePhiGradient
const ADTemplateVariableGradient< OutputType > & adGradSlnNeighbor() const override
AD grad neighbor solution getter.
bool isNodal() const override
Is this variable nodal.
const FieldVariableGradient & gradSlnNeighbor() const override
neighbor solution gradients
void evaluateOnElement(const ElemQpArg &elem_qp, const StateArg &state, bool cache_eligible) const
Evaluate solution and gradient for the elem_qp argument.
void computeIncrementAtNode(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change at the current node based on increment_vec.
virtual void computeNeighborValuesFace() override
Compute values at facial quadrature points for the neighbor.
const Elem *const & neighbor() const
Current neighboring element.
typename MooseVariableDataBase< OutputType >::ADDofValue ADDofValue
Eigen::Matrix< Real, Eigen::Dynamic, Moose::dim > RealVectorArrayValue
Definition: MooseTypes.h:149
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:343
const MooseArray< OutputType > & nodalValueArray() const override
Methods for retrieving values of variables at the nodes in a MooseArray for AuxKernelBase.
const FieldVariableValue & slnOlder() const override
const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s shape functions on a neighboring ele...
virtual std::size_t phiSize() const final
Return phi size.
const dof_id_type & nodalDofIndexNeighbor() const override
const OutputType & nodalValueDotOldNeighbor() const
virtual void insert(libMesh::NumericVector< libMesh::Number > &vector) override
Set the current local DOF values to the input vector.
Base class for time integrators.
const OutputType & nodalValueOldNeighbor() const
const DofValues & dofValuesPreviousNL() const override
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
const ADTemplateVariableCurl< OutputType > & adCurlSlnNeighbor() const override
AD curl neighbor solution getter.
const DofValues & dofValues() const override
dof values getters
const dof_id_type & nodalDofIndex() const override
const FieldVariableGradient & gradSlnNeighborDot() const
neighbor grad dots
const FieldVariableGradient & gradSlnDotDot() const
const ADTemplateVariableGradient< OutputType > & adGradSln() const override
AD grad solution getter.
const FieldVariableCurl & curlSlnOlder() const
const FieldVariablePhiCurl & curlPhi() const override final
Curl of the shape functions.
static InputParameters validParams()
virtual void computeNeighborValues() override
Compute values at quadrature points for the neighbor.
typename MooseVariableDataBase< OutputType >::ADDofValues ADDofValues
virtual void setLowerDofValues(const DenseVector< DofValue > &values) override
Set local DOF values for a lower dimensional element and evaluate the values on quadrature points...
const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
typename MooseVariableField< Real >::OutputShape OutputShape
std::vector< GradientType > _current_elem_qp_functor_gradient
The values of the gradient for the _current_elem_qp_functor_elem.
const FieldVariableValue & matrixTagValue(TagID tag) const override
const FieldVariablePhiDivergence & divPhiFaceNeighbor() const
bool computingSecond() const override final
Whether or not this variable is computing any second derivatives.
DofValue getElementalValueOld(const Elem *elem, unsigned int idx=0) const
Get the old value of this variable on an element.
std::vector< GradientType > _current_elem_side_qp_functor_gradient
The values of the gradient for the _current_elem_side_qp_functor_elem_side.
const FieldVariableDivergence & divSlnOldNeighbor() const
DofValue getNodalValue(const Node &node) const
Get the value of this variable at given node.
virtual const FieldVariablePhiValue & phiLower() const override
Return the variable&#39;s shape functions on a lower-dimensional element.
const FieldVariablePhiGradient & gradPhiLower() const
bool usesGradPhi() const
Whether or not this variable is actually using the shape function gradient.
const ADTemplateVariableTestGradient< OutputShape > & adGradPhi() const
const FieldVariableValue & uDotOldNeighbor() const
const VariableValue & duDotDu() const
bool supportsElemSideQpArg() const override final
Whether this functor supports evaluation with ElemSideQpArg.
const Elem * _current_elem_qp_functor_elem
Keep track of the current elem-qp functor element in order to enable local caching (e...
void evaluateOnElementSide(const ElemSideQpArg &elem_side_qp, const StateArg &state, bool cache_eligible) const
Evaluate solution and gradient for the elem_side_qp argument.
const DofValues & vectorTagDofValue(TagID tag) const override
State argument for evaluating functors.
const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
bool usesPhi() const
Whether or not this variable is actually using the shape function value.
const DofValues & dofValuesDotOld() const override
const FieldVariableGradient & gradSlnOlder() const
typename MooseVariableDataBase< OutputType >::DofValues DofValues
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
void prepare() override
Prepare the elemental degrees of freedom.
const FieldVariableDivergence & divSlnOlderNeighbor() const
const DofValues & dofValuesDotNeighbor() const override
const FieldVariableGradient & gradSlnDot() const
element gradient dots
std::vector< ValueType > _current_elem_side_qp_functor_sln
The values of the solution for the _current_elem_side_qp_functor_elem_side.
const ADDofValues & adDofValues() const override
Return the AD dof values.
const VariableValue & duDotDotDuNeighbor() const
void reinitNode() override
const Moose::ADType< OutputType >::type & adNodalValue() const
libMesh::TensorTools::DecrementRank< OutputType >::type OutputDivergence
const FieldVariablePhiGradient & gradPhi() const override final
Return the gradients of the variable&#39;s elemental shape functions.
const FieldVariableCurl & curlSlnOlderNeighbor() const
typename MooseVariableField< Real >::FieldVariableTestCurl FieldVariableTestCurl
const OutputType & nodalValueDotNeighborResidual() const
unsigned int oldestSolutionStateRequested() const override final
The oldest solution state that is requested for this variable (0 = current, 1 = old, 2 = older, etc).
const ADTemplateVariableValue< OutputType > & adSlnLower() const
lower-d element solution
typename MooseVariableField< Real >::OutputDivergence OutputDivergence
const FieldVariableDivergence & divSln() const
element divergence
std::unique_ptr< MooseVariableData< OutputType > > _element_data
Holder for all the data associated with the "main" element.
const FieldVariableValue & slnOld() const override
virtual void residualSetup() override
typename MooseVariableDataBase< OutputType >::DofValue DofValue
typename MooseVariableField< Real >::FieldVariablePhiDivergence FieldVariablePhiDivergence
const DofValues & dofValuesPreviousNLNeighbor() const override
DofValue getElementalValue(const Elem *elem, unsigned int idx=0) const
Get the current value of this variable on an element.
const Elem *const & currentElem() const override
Current element this variable is evaluated at.
const OutputType & nodalValueOlderNeighbor() const
const DofValues & dofValuesOldNeighbor() const override
const OutputType & nodalValueDotOld() const
const FieldVariableValue & slnOldNeighbor() const override
DotType evaluateDot(const ElemQpArg &elem_qp, const StateArg &state) const override final
const DofValues & dofValuesDotDotNeighbor() const override
virtual std::size_t phiNeighborSize() const final
Return phiNeighbor size.
Argument for requesting functor evaluation at quadrature point locations on an element side...
const FieldVariableValue & slnPreviousNLNeighbor() const
virtual void computeElemValues() override
Actually compute variable values from the solution vectors.
Moose::ShapeType< OutputType >::type OutputShape
const MappedArrayVariablePhiGradient & arrayGradPhiFace() const
typename Moose::ADType< Real >::type FunctorArg
const FieldVariableDivergence & divSlnNeighbor() const
neighbor solution divergence
const FieldVariablePhiGradient & gradPhiFace() const override final
Return the gradients of the variable&#39;s shape functions on an element face.
uint8_t dof_id_type
const FieldVariablePhiDivergence & divPhiFace() const
const FieldVariablePhiValue & phiFaceNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element face.
const OutputType & nodalValuePreviousNL() const
const ADTemplateVariableValue< OutputType > & adUDotDotNeighbor() const override
AD neighbor second time derivative getter.
void prepareAux() override
void clearAllDofIndices() final