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 
72 
89 
91 
93 
94  virtual bool isFV() const override { return true; }
95 
96  // TODO: many of these functions are not relevant to FV variables but are
97  // still called at various points from existing moose code paths. Ideally we
98  // would figure out how to remove calls to these functions and then allow
99  // throwing mooseError's from them instead of silently doing nothing (e.g.
100  // reinitNodes, reinitAux, prepareLowerD, etc.).
101 
102  virtual void prepare() override final {}
103  virtual void prepareNeighbor() override final {}
104  virtual void prepareAux() override final;
105  virtual void reinitNode() override final {}
106  virtual void reinitNodes(const std::vector<dof_id_type> & /*nodes*/) override final {}
107  virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & /*nodes*/) override final {}
108  virtual void reinitAux() override final {}
109  virtual void reinitAuxNeighbor() override final {}
110  virtual void prepareLowerD() override final {}
111 
112  virtual const dof_id_type & nodalDofIndex() const override final
113  {
114  mooseError("nodalDofIndex not supported by MooseVariableFVBase");
115  }
116  virtual const dof_id_type & nodalDofIndexNeighbor() const override final
117  {
118  mooseError("nodalDofIndexNeighbor not supported by MooseVariableFVBase");
119  }
120  virtual std::size_t phiSize() const override final { return _phi.size(); }
121  virtual std::size_t phiFaceSize() const override final { return _phi_face.size(); }
122  virtual std::size_t phiNeighborSize() const override final { return _phi_neighbor.size(); }
123  virtual std::size_t phiFaceNeighborSize() const override final
124  {
125  return _phi_face_neighbor.size();
126  }
127  virtual std::size_t phiLowerSize() const override final
128  {
129  mooseError("phiLowerSize not supported by MooseVariableFVBase");
130  }
131 
132  virtual void computeElemValuesFace() override;
133  virtual void computeNeighborValuesFace() override;
134  virtual void computeNeighborValues() override;
135  virtual void computeLowerDValues() override final
136  {
137  // mooseError("computeLowerDValues not supported by MooseVariableFVBase");
138  }
139  virtual void computeNodalNeighborValues() override final
140  {
141  // mooseError("computeNodalNeighborValues not supported by MooseVariableFVBase");
142  }
143  virtual void computeNodalValues() override final
144  {
145  // mooseError("computeNodalValues not supported by MooseVariableFVBase");
146  }
147  virtual const std::vector<dof_id_type> & dofIndicesLower() const override final;
148 
149  unsigned int numberOfDofs() const override final { return _element_data->numberOfDofs(); }
150 
151  virtual unsigned int numberOfDofsNeighbor() override final
152  {
153  mooseError("numberOfDofsNeighbor not supported by MooseVariableFVBase");
154  }
155 
156  virtual bool isNodal() const override final { return false; }
157 
158  bool hasDoFsOnNodes() const override final { return false; }
159 
160  libMesh::FEContinuity getContinuity() const override final
161  {
162  return _element_data->getContinuity();
163  };
164 
165  virtual bool isNodalDefined() const override final { return false; }
166 
167  virtual void setNodalValue(const OutputType & value, unsigned int idx = 0) override;
168 
169  virtual void setDofValue(const OutputData & value, unsigned int index) override;
170 
171  void clearDofIndices() override;
172 
173  virtual void prepareIC() override;
174 
175  virtual const Elem * const & currentElem() const override { return _element_data->currentElem(); }
176 
177  virtual void getDofIndices(const Elem * elem,
178  std::vector<dof_id_type> & dof_indices) const override
179  {
180  return _element_data->getDofIndices(elem, dof_indices);
181  }
182 
183  virtual const std::vector<dof_id_type> & dofIndices() const final
184  {
185  return _element_data->dofIndices();
186  }
187  virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const final
188  {
189  return _neighbor_data->dofIndices();
190  }
191 
193 
194  void clearAllDofIndices() final;
195 
196  const DoFValue & nodalVectorTagValue(TagID) const override
197  {
198  mooseError("nodalVectorTagValue not implemented for finite volume variables.");
199  }
200  const DoFValue & nodalMatrixTagValue(TagID) const override
201  {
202  mooseError("nodalMatrixTagValue not implemented for finite volume variables.");
203  }
204 
205  const FieldVariableValue & vectorTagValue(TagID tag) const override
206  {
207  return _element_data->vectorTagValue(tag);
208  }
209  const DoFValue & vectorTagDofValue(TagID tag) const override
210  {
211  return _element_data->vectorTagDofValue(tag);
212  }
213  const FieldVariableValue & matrixTagValue(TagID tag) const override
214  {
215  return _element_data->matrixTagValue(tag);
216  }
218  {
219  return _neighbor_data->vectorTagValue(tag);
220  }
222  {
223  return _neighbor_data->matrixTagValue(tag);
224  }
225 
226  const FieldVariableValue & uDot() const { return _element_data->uDot(); }
227  const FieldVariableValue & sln() const override { return _element_data->sln(Moose::Current); }
228  const FieldVariableValue & slnOld() const override { return _element_data->sln(Moose::Old); }
229  const FieldVariableValue & slnOlder() const override { return _element_data->sln(Moose::Older); }
230  const FieldVariableGradient & gradSln() const override
231  {
232  return _element_data->gradSln(Moose::Current);
233  }
234  const FieldVariableGradient & gradSlnOld() const override
235  {
236  return _element_data->gradSln(Moose::Old);
237  }
238  const FieldVariableValue & uDotNeighbor() const { return _neighbor_data->uDot(); }
239  const FieldVariableValue & slnNeighbor() const override
240  {
241  return _neighbor_data->sln(Moose::Current);
242  }
243  const FieldVariableValue & slnOldNeighbor() const override
244  {
245  return _neighbor_data->sln(Moose::Old);
246  }
247  const FieldVariableGradient & gradSlnNeighbor() const override
248  {
249  return _neighbor_data->gradSln(Moose::Current);
250  }
251  const FieldVariableGradient & gradSlnOldNeighbor() const override
252  {
253  return _neighbor_data->gradSln(Moose::Old);
254  }
255 
256  const VariableValue & duDotDu() const { return _element_data->duDotDu(); }
257  const VariableValue & duDotDotDu() const { return _element_data->duDotDotDu(); }
258  const VariableValue & duDotDuNeighbor() const { return _neighbor_data->duDotDu(); }
259  const VariableValue & duDotDotDuNeighbor() const { return _neighbor_data->duDotDotDu(); }
260 
263  {
264  return _element_data->adSln();
265  }
267  {
268  return _element_data->adGradSln();
269  }
270 
281  virtual const VectorValue<ADReal> & adGradSln(const Elem * const elem,
282  const StateArg & state,
283  const bool correct_skewness = false) const;
284 
295  virtual VectorValue<ADReal>
296  adGradSln(const FaceInfo & fi, const StateArg & state, const bool correct_skewness = false) const;
297 
310  virtual VectorValue<ADReal> uncorrectedAdGradSln(const FaceInfo & fi,
311  const StateArg & state,
312  const bool correct_skewness = false) const;
313 
321  const StateArg & state,
322  bool correct_skewness = false) const;
323 
325  {
326  return _element_data->adSecondSln();
327  }
329  {
330  return _element_data->adUDot();
331  }
333  {
334  return _element_data->adUDotDot();
335  }
337  {
338  return _element_data->adGradSlnDot();
339  }
341  {
342  mooseError("We don't currently implement curl for FV");
343  }
344 
347  {
348  return _neighbor_data->adSln();
349  }
351  {
352  return _neighbor_data->adGradSln();
353  }
355  {
356  return _neighbor_data->adSecondSln();
357  }
359  {
360  return _neighbor_data->adUDot();
361  }
363  {
364  return _neighbor_data->adUDotDot();
365  }
367  {
368  return _neighbor_data->adGradSlnDot();
369  }
371  {
372  mooseError("We don't currently implement curl for FV");
373  }
374 
380  virtual void computeElemValues() override;
388  virtual void computeFaceValues(const FaceInfo & fi) override;
389 
393  virtual void setDofValues(const DenseVector<OutputData> & values) override;
394  virtual void setLowerDofValues(const DenseVector<OutputData> & values) override;
395 
400  OutputData getElementalValue(const Elem * elem, unsigned int idx = 0) const;
405  OutputData getElementalValueOld(const Elem * elem, unsigned int idx = 0) const;
410  OutputData getElementalValueOlder(const Elem * elem, unsigned int idx = 0) const;
411 
412  virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) override;
413  virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) override;
414  virtual void add(libMesh::NumericVector<libMesh::Number> & vector) override;
415 
416  const DoFValue & dofValues() const override;
417  const DoFValue & dofValuesOld() const override;
418  const DoFValue & dofValuesOlder() const override;
419  const DoFValue & dofValuesPreviousNL() const override;
420  const DoFValue & dofValuesNeighbor() const override;
421  const DoFValue & dofValuesOldNeighbor() const override;
422  const DoFValue & dofValuesOlderNeighbor() const override;
423  const DoFValue & dofValuesPreviousNLNeighbor() const override;
424  const DoFValue & dofValuesDot() const override;
425  const DoFValue & dofValuesDotNeighbor() const override;
426  const DoFValue & dofValuesDotOld() const override;
427  const DoFValue & dofValuesDotOldNeighbor() const override;
428  const DoFValue & dofValuesDotDot() const override;
429  const DoFValue & dofValuesDotDotNeighbor() const override;
430  const DoFValue & dofValuesDotDotOld() const override;
431  const DoFValue & dofValuesDotDotOldNeighbor() const override;
432  const MooseArray<libMesh::Number> & dofValuesDuDotDu() const override;
433  const MooseArray<libMesh::Number> & dofValuesDuDotDuNeighbor() const override;
434  const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const override;
436 
437  const MooseArray<ADReal> & adDofValues() const override;
438  const MooseArray<ADReal> & adDofValuesNeighbor() const override;
439  const MooseArray<ADReal> & adDofValuesDot() const override;
440 
443  OutputType getValue(const Elem * elem) const;
444 
449  typename OutputTools<OutputType>::OutputGradient getGradient(const Elem * elem) const;
450 
455  bool hasDirichletBC() const
456  {
457  return _element_data->hasDirichletBC() || _neighbor_data->hasDirichletBC();
458  }
459 
460  std::pair<bool, const FVDirichletBCBase *> getDirichletBC(const FaceInfo & fi) const;
461 
462  std::pair<bool, std::vector<const FVFluxBC *>> getFluxBCs(const FaceInfo & fi) const;
463 
464  virtual void residualSetup() override;
465  virtual void initialSetup() override;
466  virtual void jacobianSetup() override;
467  virtual void timestepSetup() override;
468  virtual void meshChanged() override;
469 
477  ADReal getElemValue(const Elem * elem, const StateArg & state) const;
478 
479  void setActiveTags(const std::set<TagID> & vtags) override;
480 
486  virtual void requireQpComputations() const {}
487 
500  virtual bool isDirichletBoundaryFace(const FaceInfo & fi,
501  const Elem * elem,
502  const Moose::StateArg & state) const;
503 
504  bool supportsFaceArg() const override final { return true; }
505  bool supportsElemSideQpArg() const override final { return true; }
506 
507 protected:
522  virtual ADReal getDirichletBoundaryFaceValue(const FaceInfo & fi,
523  const Elem * elem,
524  const Moose::StateArg & state) const;
525 
540  bool isExtrapolatedBoundaryFace(const FaceInfo & fi,
541  const Elem * elem,
542  const Moose::StateArg & state) const override;
543 
544 private:
548 
549  ValueType evaluate(const ElemArg & elem, const StateArg &) const override final;
550  ValueType evaluate(const FaceArg & face, const StateArg &) const override final;
551  ValueType evaluate(const NodeArg & node, const StateArg &) const override final;
552  ValueType evaluate(const ElemPointArg & elem_point, const StateArg & state) const override final;
553  ValueType evaluate(const ElemQpArg & elem_qp, const StateArg & state) const override final;
554  ValueType evaluate(const ElemSideQpArg & elem_side_qp,
555  const StateArg & state) const override final;
556  GradientType evaluateGradient(const ElemQpArg & qp_arg, const StateArg &) const override final;
557  GradientType evaluateGradient(const ElemArg & elem_arg, const StateArg &) const override final;
558  GradientType evaluateGradient(const FaceArg & face, const StateArg &) const override final;
559  DotType evaluateDot(const ElemArg & elem, const StateArg &) const override final;
560  DotType evaluateDot(const FaceArg & face, const StateArg &) const override final;
561  DotType evaluateDot(const ElemQpArg & elem_qp, const StateArg &) const override final;
562 
567 
569  bool _dirichlet_map_setup = false;
570 
575 
577  bool _flux_map_setup = false;
578 
579 public:
580  const MooseArray<OutputType> & nodalValueArray() const override
581  {
582  mooseError("Finite volume variables do not have defined values at nodes.");
583  }
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  }
592 
593  bool computingSecond() const override final { return false; }
594  bool computingCurl() const override final { return false; }
595  bool computingDiv() const override final { return false; }
596  bool usesSecondPhiNeighbor() const override final { return false; }
597 
598  const FieldVariablePhiValue & phi() const override final { return _phi; }
599  const FieldVariablePhiGradient & gradPhi() const override final { return _grad_phi; }
600  const FieldVariablePhiSecond & secondPhi() const override final
601  {
602  mooseError("We don't currently implement second derivatives for FV");
603  }
604  const FieldVariablePhiValue & curlPhi() const override final
605  {
606  mooseError("We don't currently implement curl for FV");
607  }
608  const FieldVariablePhiDivergence & divPhi() const override final
609  {
610  mooseError("We don't currently implement divergence for FV");
611  }
612 
613  const FieldVariablePhiValue & phiFace() const override final { return _phi_face; }
614  const FieldVariablePhiGradient & gradPhiFace() const override final { return _grad_phi_face; }
615  const FieldVariablePhiSecond & secondPhiFace() const override final
616  {
617  mooseError("We don't currently implement second derivatives for FV");
618  }
619 
620  const FieldVariablePhiValue & phiFaceNeighbor() const override final
621  {
622  return _phi_face_neighbor;
623  }
624  const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
625  {
627  }
628  const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final
629  {
630  mooseError("We don't currently implement second derivatives for FV");
631  }
632 
633  const FieldVariablePhiValue & phiNeighbor() const override final { return _phi_neighbor; }
634  const FieldVariablePhiGradient & gradPhiNeighbor() const override final
635  {
636  return _grad_phi_neighbor;
637  }
638  const FieldVariablePhiSecond & secondPhiNeighbor() const override final
639  {
640  mooseError("We don't currently implement second derivatives for FV");
641  }
642 
643  virtual const FieldVariablePhiValue & phiLower() const override;
644 
645  unsigned int oldestSolutionStateRequested() const override final;
646 
669  bool two_term_expansion,
670  bool correct_skewness,
671  const Elem * elem_side_to_extrapolate_from,
672  const StateArg & state) const;
673 
676 
677 protected:
681  void clearCaches();
682 
684 
686  std::unique_ptr<MooseVariableDataFV<OutputType>> _element_data;
687 
689  std::unique_ptr<MooseVariableDataFV<OutputType>> _neighbor_data;
690 
691 private:
696 
706 
709  mutable const Elem * _prev_elem;
710 
713  std::unordered_map<BoundaryID, const FVDirichletBCBase *> _boundary_id_to_dirichlet_bc;
714 
717  std::unordered_map<BoundaryID, std::vector<const FVFluxBC *>> _boundary_id_to_flux_bc;
718 
722  [[noreturn]] void lowerDError() const;
723 
724 protected:
726  mutable std::unordered_map<const Elem *, VectorValue<ADReal>> _elem_to_grad;
727 
730 
733  mutable VectorValue<ADReal> _temp_cell_gradient;
734 
737 
742 
743  friend void Moose::initDofIndices<>(MooseVariableFV<OutputType> &, const Elem &);
744 };
745 
746 template <typename OutputType>
747 inline const MooseArray<ADReal> &
749 {
750  return _element_data->adDofValues();
751 }
752 
753 template <typename OutputType>
754 inline const MooseArray<ADReal> &
756 {
757  return _neighbor_data->adDofValues();
758 }
759 
760 template <typename OutputType>
761 inline const MooseArray<ADReal> &
763 {
764  return _element_data->adDofValuesDot();
765 }
766 
767 template <typename OutputType>
769 MooseVariableFV<OutputType>::evaluate(const ElemArg & elem_arg, const StateArg & state) const
770 {
771  return getElemValue(elem_arg.elem, state);
772 }
773 
774 template <typename OutputType>
776 MooseVariableFV<OutputType>::evaluate(const ElemPointArg & elem_point, const StateArg & state) const
777 {
778  return (*this)(elem_point.makeElem(), state) +
779  (elem_point.point - elem_point.elem->vertex_average()) *
780  this->gradient(elem_point.makeElem(), state);
781 }
782 
783 template <typename OutputType>
785 MooseVariableFV<OutputType>::evaluate(const ElemQpArg & elem_qp, const StateArg & state) const
786 {
787  return (*this)(ElemPointArg{elem_qp.elem,
788  elem_qp.point,
789  _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage},
790  state);
791 }
792 
793 template <typename OutputType>
796  const StateArg & state) const
797 {
798  return (*this)(ElemPointArg{elem_side_qp.elem,
799  elem_side_qp.point,
800  _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage},
801  state);
802 }
803 
804 template <typename OutputType>
807  const StateArg & state) const
808 {
809  return adGradSln(
810  qp_arg.elem, state, _face_interp_method == Moose::FV::InterpMethod::SkewCorrectedAverage);
811 }
812 
813 template <typename OutputType>
816  const StateArg & state) const
817 {
818  return adGradSln(elem_arg.elem, state, elem_arg.correct_skewness);
819 }
820 
821 template <typename OutputType>
824 {
825  mooseAssert(face.fi, "We must have a non-null face information");
826  return adGradSln(*face.fi, state, face.correct_skewness);
827 }
828 
829 template <typename OutputType>
830 void
831 MooseVariableFV<OutputType>::setActiveTags(const std::set<TagID> & vtags)
832 {
833  _element_data->setActiveTags(vtags);
834  _neighbor_data->setActiveTags(vtags);
835 }
836 
837 template <typename OutputType>
838 const std::vector<dof_id_type> &
840 {
841  static const std::vector<dof_id_type> empty;
842  return empty;
843 }
844 
845 template <typename OutputType>
846 void
848 {
849  determineBoundaryToDirichletBCMap();
850  determineBoundaryToFluxBCMap();
852 }
853 
854 template <typename OutputType>
855 void
857 {
858  _prev_elem = nullptr;
859  _dirichlet_map_setup = false;
860  _flux_map_setup = false;
862 }
863 
864 template <typename OutputType>
865 void
867 {
868  _dirichlet_map_setup = false;
869  _flux_map_setup = false;
871 }
872 
873 template <typename OutputType>
876 {
877  lowerDError();
878 }
879 
880 template <typename OutputType>
881 void
883 {
884  mooseError("Lower dimensional element support not implemented for finite volume variables");
885 }
886 
887 // Declare all the specializations, as the template specialization declaration below must know
888 template <>
889 ADReal MooseVariableFV<Real>::evaluateDot(const ElemArg & elem, const StateArg & state) const;
890 template <>
891 ADReal MooseVariableFV<Real>::evaluateDot(const FaceArg & elem_arg, const StateArg & state) const;
892 template <>
893 ADReal MooseVariableFV<Real>::evaluateDot(const ElemQpArg & elem_arg, const StateArg & state) const;
894 
895 // Prevent implicit instantiation in other translation units where these classes are used
896 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:609
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 DoFValue & dofValues() const override
dof values getters
const ADTemplateVariableSecond< OutputType > & adSecondSln() const override
AD second solution getter.
typename OutputTools< typename Moose::ADType< T >::type >::VariableCurl ADTemplateVariableCurl
Definition: MooseTypes.h:611
const MooseArray< ADReal > & adDofValuesNeighbor() const override
Return the AD neighbor dof values.
Moose::FaceArg FaceArg
const bool _cache_cell_gradients
Whether to cache cell gradients.
Moose::NodeArg NodeArg
const DoFValue & dofValuesPreviousNLNeighbor() const override
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.
virtual void computeNeighborValues() override
Compute values at quadrature points for the neighbor.
const DoFValue & dofValuesOld() const override
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:210
Class for stuff related to variables.
virtual void meshChanged()
Called on this object when the mesh changes.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
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.
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
const MooseArray< ADReal > & adDofValues() const override
Return the AD dof values.
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:603
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
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 DoFValue & dofValuesDotDotOld() 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...
const DoFValue & dofValuesDotOld() const override
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
A structure for storing the various lists that contain the names of the items to be exported...
void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
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.
Moose::DOFType< Real >::type OutputData
const DoFValue & dofValuesDotDot() const override
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.
const DoFValue & nodalMatrixTagValue(TagID) const override
OutputData getElementalValueOld(const Elem *elem, unsigned int idx=0) const
Get the old value of this variable on an element.
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
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.
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
const DoFValue & dofValuesDotNeighbor() const override
const FieldVariablePhiValue & phi() const override final
Return the variable&#39;s elemental shape functions.
Point point
The physical location of the quadrature point.
const DoFValue & dofValuesOlderNeighbor() 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 jacobianSetup() override
const MooseArray< libMesh::Number > & dofValuesDuDotDotDuNeighbor() const override
virtual void setLowerDofValues(const DenseVector< OutputData > &values) override
Set local DOF values for a lower dimensional element and evaluate the values on quadrature points...
const ADTemplateVariableSecond< OutputType > & adSecondSlnNeighbor() const override
AD second neighbor solution getter.
Moose::ElemPointArg ElemPointArg
const DoFValue & dofValuesPreviousNL() const override
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:606
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.
typename MooseVariableField< OutputType >::DoFValue DoFValue
const DoFValue & dofValuesOldNeighbor() const override
Argument for requesting functor evaluation at a quadrature point location in an element.
const libMesh::Elem * elem
The element.
virtual void setDofValue(const OutputData &value, unsigned int index) override
Degree of freedom value setters.
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
virtual std::size_t phiLowerSize() const override final
Return the number of shape functions on the lower dimensional element for this variable.
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:314
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 MooseArray< ADReal > & adDofValuesDot() const override
Return the AD time derivatives at dofs.
const FieldVariableGradient & gradSln() const override
element gradients
void lowerDError() const
Emit an error message for unsupported lower-d ops.
const DoFValue & dofValuesDotOldNeighbor() const override
const VariableValue & duDotDu() const
const DoFValue & dofValuesOlder() const override
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 DoFValue & dofValuesDotDotOldNeighbor() const override
const VariableValue & duDotDotDuNeighbor() const
virtual void setNodalValue(const OutputType &value, unsigned int idx=0) override
const FieldVariableValue & uDotNeighbor() const
const DoFValue & nodalVectorTagValue(TagID) const override
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.
const DoFValue & vectorTagDofValue(TagID tag) const override
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 DoFValue & dofValuesNeighbor() 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
static InputParameters validParams()
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
bool supportsElemSideQpArg() const override final
Whether this functor supports evaluation with ElemSideQpArg.
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
const InputParameters & parameters() const
Get the parameters of the object.
State argument for evaluating functors.
const FieldVariableValue & vectorTagValue(TagID tag) const override
tag values getters
const FieldVariablePhiGradient & _grad_phi_neighbor
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...
OutputData getElementalValue(const Elem *elem, unsigned int idx=0) const
Get the current value of this variable on an element.
const FieldVariablePhiValue & _phi
Shape functions.
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
virtual void setDofValues(const DenseVector< OutputData > &values) override
Set local DOF values and evaluate the values on quadrature points.
const MooseArray< OutputType > & nodalValueOlderArray() const override
virtual void timestepSetup() 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
OutputData getElementalValueOlder(const Elem *elem, unsigned int idx=0) const
Get the older value of this variable on an element.
const DoFValue & dofValuesDot() const override
Moose::ElemQpArg ElemQpArg
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...
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
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
const DoFValue & dofValuesDotDotNeighbor() const override
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.