https://mooseframework.inl.gov
MooseVariableFV.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 "MooseVariableField.h"
14 #include "SubProblem.h"
15 #include "MooseMesh.h"
16 #include "MooseVariableDataFV.h"
17 
18 #include "libmesh/numeric_vector.h"
19 #include "libmesh/dof_map.h"
20 #include "libmesh/elem.h"
21 #include "libmesh/quadrature.h"
22 #include "libmesh/dense_vector.h"
23 #include "libmesh/enum_fe_family.h"
24 
25 template <typename>
27 
29 class FVDirichletBCBase;
30 class FVFluxBC;
31 
32 namespace libMesh
33 {
34 template <typename>
35 class NumericVector;
36 }
37 
51 template <typename OutputType>
52 class MooseVariableFV : public MooseVariableField<OutputType>
53 {
54 public:
58 
64 
69 
74 
91 
93 
95 
96  virtual bool isFV() const override { return true; }
97 
98  // TODO: many of these functions are not relevant to FV variables but are
99  // still called at various points from existing moose code paths. Ideally we
100  // would figure out how to remove calls to these functions and then allow
101  // throwing mooseError's from them instead of silently doing nothing (e.g.
102  // reinitNodes, reinitAux, prepareLowerD, etc.).
103 
104  virtual void prepare() override final {}
105  virtual void prepareNeighbor() override final {}
106  virtual void prepareAux() override final;
107  virtual void reinitNode() override final {}
108  virtual void reinitNodes(const std::vector<dof_id_type> & /*nodes*/) override final {}
109  virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & /*nodes*/) override final {}
110  virtual void reinitAux() override final {}
111  virtual void reinitAuxNeighbor() override final {}
112  virtual void prepareLowerD() override final {}
113 
114  virtual const dof_id_type & nodalDofIndex() const override final
115  {
116  mooseError("nodalDofIndex not supported by MooseVariableFVBase");
117  }
118  virtual const dof_id_type & nodalDofIndexNeighbor() const override final
119  {
120  mooseError("nodalDofIndexNeighbor not supported by MooseVariableFVBase");
121  }
122  virtual std::size_t phiSize() const override final { return _phi.size(); }
123  virtual std::size_t phiFaceSize() const override final { return _phi_face.size(); }
124  virtual std::size_t phiNeighborSize() const override final { return _phi_neighbor.size(); }
125  virtual std::size_t phiFaceNeighborSize() const override final
126  {
127  return _phi_face_neighbor.size();
128  }
129  virtual std::size_t phiLowerSize() const override final
130  {
131  mooseError("phiLowerSize not supported by MooseVariableFVBase");
132  }
133 
134  virtual void computeElemValuesFace() override;
135  virtual void computeNeighborValuesFace() override;
136  virtual void computeNeighborValues() override;
137  virtual void computeLowerDValues() override final
138  {
139  // mooseError("computeLowerDValues not supported by MooseVariableFVBase");
140  }
141  virtual void computeNodalNeighborValues() override final
142  {
143  // mooseError("computeNodalNeighborValues not supported by MooseVariableFVBase");
144  }
145  virtual void computeNodalValues() override final
146  {
147  // mooseError("computeNodalValues not supported by MooseVariableFVBase");
148  }
149  virtual const std::vector<dof_id_type> & dofIndicesLower() const override final;
150 
151  unsigned int numberOfDofs() const override final { return _element_data->numberOfDofs(); }
152 
153  virtual unsigned int numberOfDofsNeighbor() override final
154  {
155  mooseError("numberOfDofsNeighbor not supported by MooseVariableFVBase");
156  }
157 
158  virtual bool isNodal() const override final { return false; }
159 
160  bool hasDoFsOnNodes() const override final { return false; }
161 
162  libMesh::FEContinuity getContinuity() const override final
163  {
164  return _element_data->getContinuity();
165  };
166 
167  virtual bool isNodalDefined() const override final { return false; }
168 
169  virtual void setNodalValue(const OutputType & value, unsigned int idx = 0) override;
170 
171  virtual void setDofValue(const DofValue & value, unsigned int index) override;
172 
173  void clearDofIndices() override;
174 
175  virtual void prepareIC() override;
176 
177  virtual const Elem * const & currentElem() const override { return _element_data->currentElem(); }
178 
179  virtual void getDofIndices(const Elem * elem,
180  std::vector<dof_id_type> & dof_indices) const override
181  {
182  return _element_data->getDofIndices(elem, dof_indices);
183  }
184 
185  virtual const std::vector<dof_id_type> & dofIndices() const final
186  {
187  return _element_data->dofIndices();
188  }
189  virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const final
190  {
191  return _neighbor_data->dofIndices();
192  }
193 
195 
196  void clearAllDofIndices() final;
197 
198  const DofValues & nodalVectorTagValue(TagID) const override
199  {
200  mooseError("nodalVectorTagValue not implemented for finite volume variables.");
201  }
202  const DofValues & nodalMatrixTagValue(TagID) const override
203  {
204  mooseError("nodalMatrixTagValue not implemented for finite volume variables.");
205  }
206 
207  const FieldVariableValue & vectorTagValue(TagID tag) const override
208  {
209  return _element_data->vectorTagValue(tag);
210  }
211  const DofValues & vectorTagDofValue(TagID tag) const override
212  {
213  return _element_data->vectorTagDofValue(tag);
214  }
215  const FieldVariableValue & matrixTagValue(TagID tag) const override
216  {
217  return _element_data->matrixTagValue(tag);
218  }
220  {
221  return _neighbor_data->vectorTagValue(tag);
222  }
224  {
225  return _neighbor_data->matrixTagValue(tag);
226  }
227 
228  const FieldVariableValue & uDot() const { return _element_data->uDot(); }
229  const FieldVariableValue & sln() const override { return _element_data->sln(Moose::Current); }
230  const FieldVariableValue & slnOld() const override { return _element_data->sln(Moose::Old); }
231  const FieldVariableValue & slnOlder() const override { return _element_data->sln(Moose::Older); }
232  const FieldVariableGradient & gradSln() const override
233  {
234  return _element_data->gradSln(Moose::Current);
235  }
236  const FieldVariableGradient & gradSlnOld() const override
237  {
238  return _element_data->gradSln(Moose::Old);
239  }
240  const FieldVariableValue & uDotNeighbor() const { return _neighbor_data->uDot(); }
241  const FieldVariableValue & slnNeighbor() const override
242  {
243  return _neighbor_data->sln(Moose::Current);
244  }
245  const FieldVariableValue & slnOldNeighbor() const override
246  {
247  return _neighbor_data->sln(Moose::Old);
248  }
249  const FieldVariableGradient & gradSlnNeighbor() const override
250  {
251  return _neighbor_data->gradSln(Moose::Current);
252  }
253  const FieldVariableGradient & gradSlnOldNeighbor() const override
254  {
255  return _neighbor_data->gradSln(Moose::Old);
256  }
257 
258  const VariableValue & duDotDu() const { return _element_data->duDotDu(); }
259  const VariableValue & duDotDotDu() const { return _element_data->duDotDotDu(); }
260  const VariableValue & duDotDuNeighbor() const { return _neighbor_data->duDotDu(); }
261  const VariableValue & duDotDotDuNeighbor() const { return _neighbor_data->duDotDotDu(); }
262 
265  {
266  return _element_data->adSln();
267  }
269  {
270  return _element_data->adGradSln();
271  }
272 
283  virtual const VectorValue<ADReal> & adGradSln(const Elem * const elem,
284  const StateArg & state,
285  const bool correct_skewness = false) const;
286 
297  virtual VectorValue<ADReal>
298  adGradSln(const FaceInfo & fi, const StateArg & state, const bool correct_skewness = false) const;
299 
312  virtual VectorValue<ADReal> uncorrectedAdGradSln(const FaceInfo & fi,
313  const StateArg & state,
314  const bool correct_skewness = false) const;
315 
323  const StateArg & state,
324  bool correct_skewness = false) const;
325 
327  {
328  return _element_data->adSecondSln();
329  }
331  {
332  return _element_data->adUDot();
333  }
335  {
336  return _element_data->adUDotDot();
337  }
339  {
340  return _element_data->adGradSlnDot();
341  }
343  {
344  mooseError("We don't currently implement curl for FV");
345  }
346 
349  {
350  return _neighbor_data->adSln();
351  }
353  {
354  return _neighbor_data->adGradSln();
355  }
357  {
358  return _neighbor_data->adSecondSln();
359  }
361  {
362  return _neighbor_data->adUDot();
363  }
365  {
366  return _neighbor_data->adUDotDot();
367  }
369  {
370  return _neighbor_data->adGradSlnDot();
371  }
373  {
374  mooseError("We don't currently implement curl for FV");
375  }
376 
382  virtual void computeElemValues() override;
390  virtual void computeFaceValues(const FaceInfo & fi) override;
391 
395  virtual void setDofValues(const DenseVector<DofValue> & values) override;
396  virtual void setLowerDofValues(const DenseVector<DofValue> & values) override;
397 
402  DofValue getElementalValue(const Elem * elem, unsigned int idx = 0) const;
407  DofValue getElementalValueOld(const Elem * elem, unsigned int idx = 0) const;
412  DofValue getElementalValueOlder(const Elem * elem, unsigned int idx = 0) const;
413 
414  virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) override;
415  virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) override;
416  virtual void add(libMesh::NumericVector<libMesh::Number> & vector) override;
417 
418  const DofValues & dofValues() const override;
419  const DofValues & dofValuesOld() const override;
420  const DofValues & dofValuesOlder() const override;
421  const DofValues & dofValuesPreviousNL() const override;
422  const DofValues & dofValuesNeighbor() const override;
423  const DofValues & dofValuesOldNeighbor() const override;
424  const DofValues & dofValuesOlderNeighbor() const override;
425  const DofValues & dofValuesPreviousNLNeighbor() const override;
426  const DofValues & dofValuesDot() const override;
427  const DofValues & dofValuesDotNeighbor() const override;
428  const DofValues & dofValuesDotOld() const override;
429  const DofValues & dofValuesDotOldNeighbor() const override;
430  const DofValues & dofValuesDotDot() const override;
431  const DofValues & dofValuesDotDotNeighbor() const override;
432  const DofValues & dofValuesDotDotOld() const override;
433  const DofValues & dofValuesDotDotOldNeighbor() const override;
434  const MooseArray<libMesh::Number> & dofValuesDuDotDu() const override;
435  const MooseArray<libMesh::Number> & dofValuesDuDotDuNeighbor() const override;
436  const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const override;
438 
439  const ADDofValues & adDofValues() const override;
440  const ADDofValues & adDofValuesNeighbor() const override;
441  const ADDofValues & adDofValuesDot() const override;
442 
445  OutputType getValue(const Elem * elem) const;
446 
451  typename OutputTools<OutputType>::OutputGradient getGradient(const Elem * elem) const;
452 
457  bool hasDirichletBC() const
458  {
459  return _element_data->hasDirichletBC() || _neighbor_data->hasDirichletBC();
460  }
461 
462  std::pair<bool, const FVDirichletBCBase *> getDirichletBC(const FaceInfo & fi) const;
463 
464  std::pair<bool, std::vector<const FVFluxBC *>> getFluxBCs(const FaceInfo & fi) const;
465 
466  virtual void residualSetup() override;
467  virtual void initialSetup() override;
468  virtual void jacobianSetup() override;
469  virtual void timestepSetup() override;
470  virtual void meshChanged() override;
471 
479  ADReal getElemValue(const Elem * elem, const StateArg & state) const;
480 
481  void setActiveTags(const std::set<TagID> & vtags) override;
482 
488  virtual void requireQpComputations() const {}
489 
502  virtual bool isDirichletBoundaryFace(const FaceInfo & fi,
503  const Elem * elem,
504  const Moose::StateArg & state) const;
505 
506  bool supportsFaceArg() const override final { return true; }
507  bool supportsElemSideQpArg() const override final { return true; }
508 
509  virtual void sizeMatrixTagData() override;
510 
511 protected:
526  virtual ADReal getDirichletBoundaryFaceValue(const FaceInfo & fi,
527  const Elem * elem,
528  const Moose::StateArg & state) const;
529 
544  bool isExtrapolatedBoundaryFace(const FaceInfo & fi,
545  const Elem * elem,
546  const Moose::StateArg & state) const override;
547 
548 private:
552 
553  ValueType evaluate(const ElemArg & elem, const StateArg &) const override final;
554  ValueType evaluate(const FaceArg & face, const StateArg &) const override final;
555  ValueType evaluate(const NodeArg & node, const StateArg &) const override final;
556  ValueType evaluate(const ElemPointArg & elem_point, const StateArg & state) const override final;
557  ValueType evaluate(const ElemQpArg & elem_qp, const StateArg & state) const override final;
558  ValueType evaluate(const ElemSideQpArg & elem_side_qp,
559  const StateArg & state) const override final;
560  GradientType evaluateGradient(const ElemQpArg & qp_arg, const StateArg &) const override final;
561  GradientType evaluateGradient(const ElemArg & elem_arg, const StateArg &) const override final;
562  GradientType evaluateGradient(const FaceArg & face, const StateArg &) const override final;
563  DotType evaluateDot(const ElemArg & elem, const StateArg &) const override final;
564  DotType evaluateDot(const FaceArg & face, const StateArg &) const override final;
565  DotType evaluateDot(const ElemQpArg & elem_qp, const StateArg &) const override final;
566 
571 
573  bool _dirichlet_map_setup = false;
574 
579 
581  bool _flux_map_setup = false;
582 
583 public:
584  const MooseArray<OutputType> & nodalValueArray() const override
585  {
586  mooseError("Finite volume variables do not have defined values at nodes.");
587  }
589  {
590  mooseError("Finite volume variables do not have defined values at nodes.");
591  }
593  {
594  mooseError("Finite volume variables do not have defined values at nodes.");
595  }
596 
597  bool computingSecond() const override final { return false; }
598  bool computingCurl() const override final { return false; }
599  bool computingDiv() const override final { return false; }
600  bool usesSecondPhiNeighbor() const override final { return false; }
601 
602  const FieldVariablePhiValue & phi() const override final { return _phi; }
603  const FieldVariablePhiGradient & gradPhi() const override final { return _grad_phi; }
604  const FieldVariablePhiSecond & secondPhi() const override final
605  {
606  mooseError("We don't currently implement second derivatives for FV");
607  }
608  const FieldVariablePhiValue & curlPhi() const override final
609  {
610  mooseError("We don't currently implement curl for FV");
611  }
612  const FieldVariablePhiDivergence & divPhi() const override final
613  {
614  mooseError("We don't currently implement divergence for FV");
615  }
616 
617  const FieldVariablePhiValue & phiFace() const override final { return _phi_face; }
618  const FieldVariablePhiGradient & gradPhiFace() const override final { return _grad_phi_face; }
619  const FieldVariablePhiSecond & secondPhiFace() const override final
620  {
621  mooseError("We don't currently implement second derivatives for FV");
622  }
623 
624  const FieldVariablePhiValue & phiFaceNeighbor() const override final
625  {
626  return _phi_face_neighbor;
627  }
628  const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
629  {
631  }
632  const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final
633  {
634  mooseError("We don't currently implement second derivatives for FV");
635  }
636 
637  const FieldVariablePhiValue & phiNeighbor() const override final { return _phi_neighbor; }
638  const FieldVariablePhiGradient & gradPhiNeighbor() const override final
639  {
640  return _grad_phi_neighbor;
641  }
642  const FieldVariablePhiSecond & secondPhiNeighbor() const override final
643  {
644  mooseError("We don't currently implement second derivatives for FV");
645  }
646 
647  virtual const FieldVariablePhiValue & phiLower() const override;
648 
649  unsigned int oldestSolutionStateRequested() const override final;
650 
673  bool two_term_expansion,
674  bool correct_skewness,
675  const Elem * elem_side_to_extrapolate_from,
676  const StateArg & state) const;
677 
680 
681 protected:
685  void clearCaches();
686 
688 
690  std::unique_ptr<MooseVariableDataFV<OutputType>> _element_data;
691 
693  std::unique_ptr<MooseVariableDataFV<OutputType>> _neighbor_data;
694 
695 private:
700 
710 
713  mutable const Elem * _prev_elem;
714 
717  std::unordered_map<BoundaryID, const FVDirichletBCBase *> _boundary_id_to_dirichlet_bc;
718 
721  std::unordered_map<BoundaryID, std::vector<const FVFluxBC *>> _boundary_id_to_flux_bc;
722 
726  [[noreturn]] void lowerDError() const;
727 
728 protected:
730  mutable std::unordered_map<const Elem *, VectorValue<ADReal>> _elem_to_grad;
731 
734 
737  mutable VectorValue<ADReal> _temp_cell_gradient;
738 
741 
746 
747  friend void Moose::initDofIndices<>(MooseVariableFV<OutputType> &, const Elem &);
748 };
749 
750 template <typename OutputType>
751 inline const typename MooseVariableFV<OutputType>::ADDofValues &
753 {
754  return _element_data->adDofValues();
755 }
756 
757 template <typename OutputType>
758 inline const typename MooseVariableFV<OutputType>::ADDofValues &
760 {
761  return _neighbor_data->adDofValues();
762 }
763 
764 template <typename OutputType>
765 inline const typename MooseVariableFV<OutputType>::ADDofValues &
767 {
768  return _element_data->adDofValuesDot();
769 }
770 
771 template <typename OutputType>
773 MooseVariableFV<OutputType>::evaluate(const ElemArg & elem_arg, const StateArg & state) const
774 {
775  return getElemValue(elem_arg.elem, state);
776 }
777 
778 template <typename OutputType>
780 MooseVariableFV<OutputType>::evaluate(const ElemPointArg & elem_point, const StateArg & state) const
781 {
782  return (*this)(elem_point.makeElem(), state) +
783  (elem_point.point - elem_point.elem->vertex_average()) *
784  this->gradient(elem_point.makeElem(), state);
785 }
786 
787 template <typename OutputType>
789 MooseVariableFV<OutputType>::evaluate(const ElemQpArg & elem_qp, const StateArg & state) const
790 {
791  return (*this)(ElemPointArg{elem_qp.elem,
792  elem_qp.point,
793  _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage},
794  state);
795 }
796 
797 template <typename OutputType>
800  const StateArg & state) const
801 {
802  return (*this)(ElemPointArg{elem_side_qp.elem,
803  elem_side_qp.point,
804  _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage},
805  state);
806 }
807 
808 template <typename OutputType>
811  const StateArg & state) const
812 {
813  return adGradSln(
814  qp_arg.elem, state, _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage);
815 }
816 
817 template <typename OutputType>
820  const StateArg & state) const
821 {
822  return adGradSln(elem_arg.elem, state, elem_arg.correct_skewness);
823 }
824 
825 template <typename OutputType>
828 {
829  mooseAssert(face.fi, "We must have a non-null face information");
830  return adGradSln(*face.fi, state, face.correct_skewness);
831 }
832 
833 template <typename OutputType>
834 void
835 MooseVariableFV<OutputType>::setActiveTags(const std::set<TagID> & vtags)
836 {
837  _element_data->setActiveTags(vtags);
838  _neighbor_data->setActiveTags(vtags);
839 }
840 
841 template <typename OutputType>
842 const std::vector<dof_id_type> &
844 {
845  static const std::vector<dof_id_type> empty;
846  return empty;
847 }
848 
849 template <typename OutputType>
850 void
852 {
853  determineBoundaryToDirichletBCMap();
854  determineBoundaryToFluxBCMap();
856 }
857 
858 template <typename OutputType>
859 void
861 {
862  _prev_elem = nullptr;
863  _dirichlet_map_setup = false;
864  _flux_map_setup = false;
866 }
867 
868 template <typename OutputType>
869 void
871 {
872  _dirichlet_map_setup = false;
873  _flux_map_setup = false;
875 }
876 
877 template <typename OutputType>
880 {
881  lowerDError();
882 }
883 
884 template <typename OutputType>
885 void
887 {
888  mooseError("Lower dimensional element support not implemented for finite volume variables");
889 }
890 
891 // Declare all the specializations, as the template specialization declaration below must know
892 template <>
893 ADReal MooseVariableFV<Real>::evaluateDot(const ElemArg & elem, const StateArg & state) const;
894 template <>
895 ADReal MooseVariableFV<Real>::evaluateDot(const FaceArg & elem_arg, const StateArg & state) const;
896 template <>
897 ADReal MooseVariableFV<Real>::evaluateDot(const ElemQpArg & elem_arg, const StateArg & state) const;
898 
899 // Prevent implicit instantiation in other translation units where these classes are used
900 extern template class MooseVariableFV<Real>;
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
bool _two_term_boundary_expansion
Whether to use a two term expansion for computing boundary face values.
typename OutputTools< typename Moose::ADType< T >::type >::VariableSecond ADTemplateVariableSecond
Definition: MooseTypes.h:654
libMesh::FEContinuity getContinuity() const override final
Return the continuity of this variable.
const ADTemplateVariableValue< OutputType > & adUDotDot() const override
AD second time derivative getter.
virtual bool isNodalDefined() const override final
Is this variable defined at nodes.
GradientType evaluateGradient(const ElemQpArg &qp_arg, const StateArg &) const override final
Moose::FV::InterpMethod _face_interp_method
Decides if an average or skewed corrected average is used for the face interpolation.
const FieldVariablePhiValue & phiNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element.
const ADTemplateVariableSecond< OutputType > & adSecondSln() const override
AD second solution getter.
typename OutputTools< typename Moose::ADType< T >::type >::VariableCurl ADTemplateVariableCurl
Definition: MooseTypes.h:656
Moose::FaceArg FaceArg
const bool _cache_cell_gradients
Whether to cache cell gradients.
DofValue getElementalValueOld(const Elem *elem, unsigned int idx=0) const
Get the old value of this variable on an element.
Moose::NodeArg NodeArg
void determineBoundaryToFluxBCMap()
Setup the boundary to Flux BC map.
virtual std::size_t phiFaceSize() const override final
Return phiFace size.
virtual bool isNodal() const override final
Is this variable nodal.
const ADDofValues & adDofValuesNeighbor() const override
Return the AD neighbor dof values.
const DofValues & nodalVectorTagValue(TagID) const override
virtual void computeNeighborValues() override
Compute values at quadrature points for the neighbor.
virtual void computeLowerDValues() override final
compute values at quadrature points on the lower dimensional element
virtual void prepareAux() override final
void clearDofIndices() override
Clear out the dof indices.
unsigned int TagID
Definition: MooseTypes.h:238
virtual void setDofValues(const DenseVector< DofValue > &values) override
Set local DOF values and evaluate the values on quadrature points.
Class for stuff related to variables.
virtual void meshChanged()
Called on this object when the mesh changes.
const ADDofValues & adDofValues() const override
Return the AD dof values.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
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 FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element face.
virtual void computeNodalValues() override final
Compute nodal values of this variable.
std::unique_ptr< MooseVariableDataFV< OutputType > > _neighbor_data
Holder for all the data associated with the neighbor element.
const bool & getTwoTermBoundaryExpansion() const
Function to get wether two term boundary expansion is used for the variable.
const MooseArray< libMesh::Number > & dofValuesDuDotDu() const override
const FieldVariablePhiValue & _phi_face_neighbor
virtual void reinitAux() override final
void determineBoundaryToDirichletBCMap()
Setup the boundary to Dirichlet BC map.
Moose::StateArg StateArg
virtual void meshChanged() override
Called on this object when the mesh changes.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
const DofValues & dofValuesDotDot() const override
bool computingCurl() const override final
Whether or not this variable is computing any curl quantities.
Base class for finite volume Dirichlet boundaray conditions.
virtual void prepareLowerD() override final
Prepare a lower dimensional element&#39;s degrees of freedom.
bool correct_skewness
Whether to perform skew correction.
unsigned int numberOfDofs() const override final
Get the number of local DoFs.
std::unordered_map< BoundaryID, std::vector< const FVFluxBC * > > _boundary_id_to_flux_bc
Map from boundary ID to flux boundary conditions.
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:648
Moose::ElemArg ElemArg
virtual std::size_t phiFaceNeighborSize() const override final
Return phiFaceNeighbor size.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const FieldVariableGradient & gradSlnOldNeighbor() const override
virtual void sizeMatrixTagData() override
Size data structures related to matrix tagging.
const FieldVariablePhiGradient & _grad_phi_face_neighbor
const ADTemplateVariableValue< OutputType > & adSln() const override
AD.
OutputTools< OutputType >::OutputGradient getGradient(const Elem *elem) const
Compute the variable gradient value at a point on an element.
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...
virtual void reinitNodes(const std::vector< dof_id_type > &) override final
ValueType evaluate(const ElemArg &elem, const StateArg &) const override final
Evaluate the functor with a given element.
const FieldVariablePhiGradient & gradPhi() const override final
Return the gradients of the variable&#39;s elemental shape functions.
const DofValues & dofValuesOldNeighbor() const override
const FieldVariablePhiSecond & secondPhi() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s elemental shape functions.
Moose::ElemSideQpArg ElemSideQpArg
virtual void prepare() override final
Prepare the elemental degrees of freedom.
typename MooseVariableField< OutputType >::FieldVariablePhiGradient FieldVariablePhiGradient
const VariableValue & duDotDotDu() const
const FieldVariablePhiGradient & gradPhiFace() const override final
Return the gradients of the variable&#39;s shape functions on an element face.
A structure that is used to evaluate Moose functors at an arbitrary physical point contained within a...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
const ADTemplateVariableGradient< OutputType > & adGradSlnNeighbor() const override
AD grad neighbor solution getter.
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...
bool usesSecondPhiNeighbor() const override final
Whether or not this variable is actually using the shape function second derivatives.
virtual void add(libMesh::NumericVector< libMesh::Number > &vector) override
Add the currently cached degree of freedom values into the provided vector.
virtual void prepareIC() override
Prepare the initial condition.
VectorValue< ADReal > _temp_cell_gradient
A member to hold the cell gradient when not caching, used to return a reference (due to expensive ADR...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:42
const ADDofValues & adDofValuesDot() const override
Return the AD time derivatives at dofs.
void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
const DofValues & dofValuesOlder() const override
const FieldVariablePhiDivergence & divPhi() const override final
Divergence of the shape functions.
const ADTemplateVariableCurl< OutputType > & adCurlSlnNeighbor() const override
AD curl neighbor solution getter.
std::pair< bool, std::vector< const FVFluxBC * > > getFluxBCs(const FaceInfo &fi) const
ADReal getElemValue(const Elem *elem, const StateArg &state) const
Get the solution value for the provided element and seed the derivative for the corresponding dof ind...
typename FunctorReturnType< Moose::ADType< OutputType >::type, FunctorEvaluationKind::Gradient >::type GradientType
This rigmarole makes it so that a user can create functors that return containers (std::vector...
Definition: MooseFunctor.h:149
const FieldVariablePhiGradient & gradPhiNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element.
virtual bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
Determine whether a specified face side is a Dirichlet boundary face.
typename MooseVariableField< OutputType >::FieldVariablePhiValue FieldVariablePhiValue
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...
const MooseArray< OutputType > & nodalValueArray() const override
Methods for retrieving values of variables at the nodes in a MooseArray for AuxKernelBase.
virtual void reinitAuxNeighbor() override final
const libMesh::NumericVector< libMesh::Number > *const & _solution
The current (ghosted) solution.
libMesh::TensorTools::DecrementRank< OutputShape >::type OutputShapeDivergence
const FieldVariableValue & slnOlder() const override
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
std::pair< bool, const FVDirichletBCBase * > getDirichletBC(const FaceInfo &fi) const
MooseVariableFV(const InputParameters &parameters)
const ADTemplateVariableGradient< OutputType > & adGradSlnDot() const override
AD grad of time derivative solution getter.
OutputType getValue(const Elem *elem) const
Note: const monomial is always the case - higher order solns are reconstructed - so this is simpler f...
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const override
const MooseArray< libMesh::Number > & dofValuesDuDotDuNeighbor() const override
const FieldVariableValue & sln() const override
void clearCaches()
clear finite volume caches
const DofValues & dofValuesDotNeighbor() const override
virtual const Elem *const & currentElem() const override
Current element this variable is evaluated at.
const FieldVariableGradient & gradSlnNeighbor() const override
neighbor solution gradients
virtual VectorValue< ADReal > uncorrectedAdGradSln(const FaceInfo &fi, const StateArg &state, const bool correct_skewness=false) const
Retrieve (or potentially compute) the uncorrected gradient on the provided face.
bool hasDirichletBC() const
Returns true if a Dirichlet BC exists on the current face.
bool _dirichlet_map_setup
Whether the boundary to Dirichlet cache map has been setup yet.
std::unique_ptr< MooseVariableDataFV< OutputType > > _element_data
Holder for all the data associated with the "main" element.
virtual const FieldVariablePhiValue & phiLower() const override
Return the variable&#39;s shape functions on a lower-dimensional element.
const DofValues & dofValuesPreviousNLNeighbor() const override
const DofValues & dofValues() const override
dof values getters
A structure defining a "face" evaluation calling argument for Moose functors.
const FieldVariableValue & matrixTagValue(TagID tag) const override
const FaceInfo * fi
a face information object which defines our location in space
const ADTemplateVariableValue< OutputType > & adSlnNeighbor() const override
neighbor AD
virtual void setDofValue(const DofValue &value, unsigned int index) override
Degree of freedom value setters.
const FieldVariablePhiValue & phi() const override final
Return the variable&#39;s elemental shape functions.
Point point
The physical location of the quadrature point.
const DofValues & dofValuesDotDotOldNeighbor() const override
const FieldVariableGradient & gradSlnOld() const override
virtual void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const override
const libMesh::Elem * elem
virtual void setLowerDofValues(const DenseVector< DofValue > &values) override
Set local DOF values for a lower dimensional element and evaluate the values on quadrature points...
virtual void jacobianSetup() override
const MooseArray< libMesh::Number > & dofValuesDuDotDotDuNeighbor() const override
const ADTemplateVariableSecond< OutputType > & adSecondSlnNeighbor() const override
AD second neighbor solution getter.
Moose::ElemPointArg ElemPointArg
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:651
const libMesh::Elem * elem
Provides an interface for computing residual contributions from finite volume numerical fluxes comput...
Definition: FVFluxBC.h:23
MooseVariableFV< Real > MooseVariableFVReal
const FieldVariableValue & matrixTagValueNeighbor(TagID tag)
const FieldVariableValue & slnNeighbor() const override
virtual void computeFaceValues(const FaceInfo &fi) override
Initializes/computes variable values from the solution vectors for the face represented by fi...
bool isExtrapolatedBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const override
Returns whether this is an extrapolated boundary face.
void clearAllDofIndices() final
A structure that is used to evaluate Moose functors logically at an element/cell center.
Argument for requesting functor evaluation at a quadrature point location in an element.
DofValue getElementalValueOlder(const Elem *elem, unsigned int idx=0) const
Get the older value of this variable on an element.
const libMesh::Elem * elem
The element.
bool supportsFaceArg() const override final
Whether this functor supports evaluation with FaceArg.
virtual std::size_t phiNeighborSize() const override final
Return phiNeighbor size.
virtual void timestepSetup() override
const ADTemplateVariableCurl< OutputType > & adCurlSln() const override
AD curl solution getter.
const FieldVariablePhiValue & _phi_neighbor
const ADTemplateVariableGradient< OutputType > & adGradSlnNeighborDot() const override
AD grad of time derivative neighbor solution getter.
virtual ADReal getExtrapolatedBoundaryFaceValue(const FaceInfo &fi, bool two_term_expansion, bool correct_skewness, const Elem *elem_side_to_extrapolate_from, const StateArg &state) const
Retrieves an extrapolated boundary value for the provided face.
const FieldVariableValue & uDot() const
const FieldVariablePhiValue & curlPhi() const override final
Curl of the shape functions.
virtual const dof_id_type & nodalDofIndex() const override final
virtual bool isFV() const override
const DofValues & vectorTagDofValue(TagID tag) const override
const DofValues & dofValuesDotDotOld() const override
virtual std::size_t phiLowerSize() const override final
Return the number of shape functions on the lower dimensional element for this variable.
typename MooseVariableDataBase< OutputType >::ADDofValue ADDofValue
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:343
const ADTemplateVariableValue< OutputType > & adUDot() const override
AD time derivative getter.
virtual std::size_t phiSize() const override final
Return phi size.
const FieldVariableValue & slnOldNeighbor() const override
const DofValues & dofValuesDot() const override
const FieldVariableGradient & gradSln() const override
element gradients
void lowerDError() const
Emit an error message for unsupported lower-d ops.
const VariableValue & duDotDu() const
bool _flux_map_setup
Whether the boundary to fluxBC cache map has been setup yet.
(gc*elem+(1-gc)*neighbor)+gradient*(rf-rf&#39;)
Moose::FV::InterpMethod faceInterpolationMethod() const
const VariableValue & duDotDotDuNeighbor() const
virtual void setNodalValue(const OutputType &value, unsigned int idx=0) override
const FieldVariableValue & uDotNeighbor() const
bool hasDoFsOnNodes() const override final
Does this variable have DoFs on nodes.
virtual void insert(libMesh::NumericVector< libMesh::Number > &vector) override
Insert the currently cached degree of freedom values into the provided vector.
typename MooseVariableDataBase< OutputType >::ADDofValues ADDofValues
const FieldVariablePhiValue & phiFaceNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element face.
const MooseArray< OutputType > & nodalValueOldArray() const override
const DofValues & dofValuesDotOldNeighbor() const override
const DofValues & dofValuesOld() const override
ADReal getBoundaryFaceValue(const FaceInfo &fi, const StateArg &state, bool correct_skewness=false) const
Retrieve the solution value at a boundary face.
const ADTemplateVariableValue< OutputType > & adUDotDotNeighbor() const override
AD neighbor second time derivative getter.
virtual void reinitNodesNeighbor(const std::vector< dof_id_type > &) override final
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
static InputParameters validParams()
bool supportsElemSideQpArg() const override final
Whether this functor supports evaluation with ElemSideQpArg.
const DofValues & dofValuesOlderNeighbor() const override
bool computingDiv() const override final
Whether or not this variable is computing any divergence quantities.
const FieldVariablePhiValue & _phi_face
virtual ADReal getDirichletBoundaryFaceValue(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
Retrieves a Dirichlet boundary value for the provided face.
virtual void reinitNode() override final
State argument for evaluating functors.
const FieldVariableValue & vectorTagValue(TagID tag) const override
tag values getters
const FieldVariablePhiGradient & _grad_phi_neighbor
typename MooseVariableDataBase< OutputType >::DofValues DofValues
virtual void computeNodalNeighborValues() override final
Compute nodal values of this variable in the neighbor.
ElemArg makeElem() const
Make a ElemArg from our data.
virtual void computeElemValues() override
Initializes/computes variable values from the solution vectors for the current element being operated...
const FieldVariablePhiValue & _phi
Shape functions.
const DofValues & nodalMatrixTagValue(TagID) const override
virtual void prepareNeighbor() override final
Prepare the neighbor element degrees of freedom.
libMesh::Point point
The physical location of the quadrature point.
libMesh::TensorTools::DecrementRank< Real >::type OutputDivergence
const FieldVariablePhiGradient & _grad_phi_face
virtual void residualSetup() override
const MooseArray< OutputType > & nodalValueOlderArray() const override
virtual void timestepSetup() override
const DofValues & dofValuesNeighbor() const override
This class provides variable solution values for other classes/objects to bind to when looping over f...
virtual void computeElemValuesFace() override
Compute values at facial quadrature points.
virtual const std::vector< dof_id_type > & dofIndicesLower() const override final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
virtual unsigned int numberOfDofsNeighbor() override final
const DofValues & dofValuesDotDotNeighbor() const override
Moose::ElemQpArg ElemQpArg
typename MooseVariableDataBase< OutputType >::DofValue DofValue
const VariableValue & duDotDuNeighbor() const
const libMesh::Elem * elem
The element.
void setActiveTags(const std::set< TagID > &vtags) override
Set the active vector tags.
const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
const Elem * _prev_elem
A member used to help determine when we can return cached data as opposed to computing new data...
const DofValues & dofValuesPreviousNL() const override
InterpMethod
This codifies a set of available ways to interpolate with elem+neighbor solution information to calcu...
Definition: MathFVUtils.h:35
virtual void computeNeighborValuesFace() override
Compute values at facial quadrature points for the neighbor.
const ADTemplateVariableValue< OutputType > & adUDotNeighbor() const override
AD neighbor time derivative getter.
const FieldVariableValue & slnOld() const override
DotType evaluateDot(const ElemArg &elem, const StateArg &) const override final
Evaluate the functor time derivative with a given element.
virtual const dof_id_type & nodalDofIndexNeighbor() const override final
unsigned int oldestSolutionStateRequested() const override final
The oldest solution state that is requested for this variable (0 = current, 1 = old, 2 = older, etc).
Argument for requesting functor evaluation at quadrature point locations on an element side...
const FieldVariablePhiGradient & _grad_phi
Moose::ShapeType< Real >::type OutputShape
const DofValues & dofValuesDotOld() const override
Point vertex_average() const
virtual void requireQpComputations() const
Request that quadrature point data be (pre)computed.
std::unordered_map< BoundaryID, const FVDirichletBCBase * > _boundary_id_to_dirichlet_bc
Map from boundary ID to Dirichlet boundary conditions.
virtual const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
uint8_t dof_id_type
DofValue getElementalValue(const Elem *elem, unsigned int idx=0) const
Get the current value of this variable on an element.
const FieldVariableValue & vectorTagValueNeighbor(TagID tag)
bool computingSecond() const override final
Whether or not this variable is computing any second derivatives.
std::unordered_map< const Elem *, VectorValue< ADReal > > _elem_to_grad
A cache for storing gradients on elements.
const ADTemplateVariableGradient< OutputType > & adGradSln() const override
AD grad solution getter.