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" 51 template <
typename OutputType>
94 virtual bool isFV()
const override {
return true; }
106 virtual void reinitNodes(
const std::vector<dof_id_type> & )
override final {}
114 mooseError(
"nodalDofIndex not supported by MooseVariableFVBase");
118 mooseError(
"nodalDofIndexNeighbor not supported by MooseVariableFVBase");
120 virtual std::size_t
phiSize() const override final {
return _phi.size(); }
129 mooseError(
"phiLowerSize not supported by MooseVariableFVBase");
147 virtual const std::vector<dof_id_type> &
dofIndicesLower() const override final;
153 mooseError(
"numberOfDofsNeighbor not supported by MooseVariableFVBase");
156 virtual bool isNodal() const override final {
return false; }
167 virtual void setNodalValue(
const OutputType & value,
unsigned int idx = 0)
override;
178 std::vector<dof_id_type> & dof_indices)
const override 183 virtual const std::vector<dof_id_type> &
dofIndices() const final
198 mooseError(
"nodalVectorTagValue not implemented for finite volume variables.");
202 mooseError(
"nodalMatrixTagValue not implemented for finite volume variables.");
281 virtual const VectorValue<ADReal> &
adGradSln(
const Elem *
const elem,
283 const bool correct_skewness =
false)
const;
295 virtual VectorValue<ADReal>
312 const bool correct_skewness =
false)
const;
322 bool correct_skewness =
false)
const;
342 mooseError(
"We don't currently implement curl for FV");
372 mooseError(
"We don't currently implement curl for FV");
393 virtual void setDofValues(
const DenseVector<OutputData> & values)
override;
443 OutputType
getValue(
const Elem * elem)
const;
555 const StateArg & state)
const override final;
582 mooseError(
"Finite volume variables do not have defined values at nodes.");
586 mooseError(
"Finite volume variables do not have defined values at nodes.");
590 mooseError(
"Finite volume variables do not have defined values at nodes.");
602 mooseError(
"We don't currently implement second derivatives for FV");
606 mooseError(
"We don't currently implement curl for FV");
610 mooseError(
"We don't currently implement divergence for FV");
617 mooseError(
"We don't currently implement second derivatives for FV");
630 mooseError(
"We don't currently implement second derivatives for FV");
640 mooseError(
"We don't currently implement second derivatives for FV");
669 bool two_term_expansion,
670 bool correct_skewness,
671 const Elem * elem_side_to_extrapolate_from,
726 mutable std::unordered_map<const Elem *, VectorValue<ADReal>>
_elem_to_grad;
746 template <
typename OutputType>
750 return _element_data->adDofValues();
753 template <
typename OutputType>
757 return _neighbor_data->adDofValues();
760 template <
typename OutputType>
764 return _element_data->adDofValuesDot();
767 template <
typename OutputType>
771 return getElemValue(elem_arg.
elem, state);
774 template <
typename OutputType>
778 return (*
this)(elem_point.
makeElem(), state) +
780 this->gradient(elem_point.
makeElem(), state);
783 template <
typename OutputType>
793 template <
typename OutputType>
804 template <
typename OutputType>
813 template <
typename OutputType>
821 template <
typename OutputType>
825 mooseAssert(face.
fi,
"We must have a non-null face information");
829 template <
typename OutputType>
833 _element_data->setActiveTags(vtags);
834 _neighbor_data->setActiveTags(vtags);
837 template <
typename OutputType>
838 const std::vector<dof_id_type> &
841 static const std::vector<dof_id_type> empty;
845 template <
typename OutputType>
849 determineBoundaryToDirichletBCMap();
850 determineBoundaryToFluxBCMap();
854 template <
typename OutputType>
858 _prev_elem =
nullptr;
859 _dirichlet_map_setup =
false;
860 _flux_map_setup =
false;
864 template <
typename OutputType>
868 _dirichlet_map_setup =
false;
869 _flux_map_setup =
false;
873 template <
typename OutputType>
880 template <
typename OutputType>
884 mooseError(
"Lower dimensional element support not implemented for finite volume variables");
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
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'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
const MooseArray< ADReal > & adDofValuesNeighbor() const override
Return the AD neighbor dof values.
const bool _cache_cell_gradients
Whether to cache cell gradients.
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.
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...
const FieldVariablePhiSecond & secondPhiFace() const override final
Return the rank-2 tensor of second derivatives of the variable's shape functions on an element face...
const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
Return the gradients of the variable'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.
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'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
virtual std::size_t phiFaceNeighborSize() const override final
Return phiFaceNeighbor size.
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'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'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'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'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
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...
const FieldVariablePhiGradient & gradPhiNeighbor() const override final
Return the gradients of the variable'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'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 ...
std::pair< bool, const FVDirichletBCBase * > getDirichletBC(const FaceInfo &fi) const
MooseVariableFV(const InputParameters ¶meters)
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'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'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
const libMesh::Elem * elem
Provides an interface for computing residual contributions from finite volume numerical fluxes comput...
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
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')
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
Moose::ADType< OutputType >::type ValueType
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'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'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...
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
usingMooseVariableFieldMembers
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.
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.