https://mooseframework.inl.gov
MooseLinearVariableFV.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"
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;
32 
33 namespace libMesh
34 {
35 template <typename>
36 class NumericVector;
37 }
38 
42 template <typename OutputType>
43 class MooseLinearVariableFV : public MooseVariableField<OutputType>
44 {
45 public:
49 
55 
60 
65 
82 
85 
86  virtual bool isFV() const override { return true; }
87 
91  virtual bool isDirichletBoundaryFace(const FaceInfo & fi) const;
92 
97 
101  virtual bool needsGradientVectorStorage() const override { return _needs_cell_gradients; }
102 
103  virtual bool isExtrapolatedBoundaryFace(const FaceInfo & fi,
104  const Elem * elem,
105  const Moose::StateArg & state) const override;
106 
111  const VectorValue<Real> gradSln(const ElemInfo & elem_info) const;
112 
119  VectorValue<Real> gradSln(const FaceInfo & fi, const StateArg & state) const;
120 
121  virtual void initialSetup() override;
122 
130  Real getElemValue(const ElemInfo & elem_info, const StateArg & state) const;
131 
137 
138  const std::unordered_map<BoundaryID, LinearFVBoundaryCondition *> & getBoundaryConditionMap()
139  {
140  return _boundary_id_to_bc;
141  }
142 
143  virtual void prepareIC() override {}
144 
145  virtual bool isNodal() const override final { return false; }
146 
147  virtual bool hasDoFsOnNodes() const override final { return false; }
148 
149  virtual bool isNodalDefined() const override final { return false; }
150 
151  virtual bool supportsFaceArg() const override final { return true; }
152  virtual bool supportsElemSideQpArg() const override final { return false; }
153 
154  virtual const Elem * const & currentElem() const override;
155 
156  virtual bool computingSecond() const override final { return false; }
157  virtual bool computingCurl() const override final { return false; }
158  virtual bool computingDiv() const override final { return false; }
159  virtual bool usesSecondPhiNeighbor() const override final { return false; }
160 
161  virtual void sizeMatrixTagData() override;
162 
163 protected:
165  [[noreturn]] void timeIntegratorError() const;
166 
168  [[noreturn]] void lowerDError() const;
169 
171  [[noreturn]] void nodalError() const;
172 
174  [[noreturn]] void adError() const;
175 
179  void cacheBoundaryBCMap();
180 
182 
185 
188 
190  const std::vector<std::unique_ptr<libMesh::NumericVector<libMesh::Number>>> & _grad_container;
191 
195  std::unique_ptr<MooseVariableDataLinearFV<OutputType>> _element_data;
196 
200  std::unique_ptr<MooseVariableDataLinearFV<OutputType>> _neighbor_data;
201 
204  std::unordered_map<BoundaryID, LinearFVBoundaryCondition *> _boundary_id_to_bc;
205 
207  const unsigned int _sys_num;
208 
209  friend void Moose::initDofIndices<>(MooseLinearVariableFV<OutputType> &, const Elem &);
210 
211 private:
215 
216  virtual ValueType evaluate(const ElemArg & elem, const StateArg &) const override final;
217  virtual ValueType evaluate(const FaceArg & face, const StateArg &) const override final;
218  virtual ValueType evaluate(const NodeArg & node, const StateArg &) const override final;
219  virtual ValueType evaluate(const ElemPointArg & elem_point,
220  const StateArg & state) const override final;
221  virtual ValueType evaluate(const ElemQpArg & elem_qp,
222  const StateArg & state) const override final;
223  virtual ValueType evaluate(const ElemSideQpArg & elem_side_qp,
224  const StateArg & state) const override final;
225  virtual GradientType evaluateGradient(const ElemQpArg & qp_arg,
226  const StateArg &) const override final;
227  virtual GradientType evaluateGradient(const ElemArg & elem_arg,
228  const StateArg &) const override final;
229  virtual GradientType evaluateGradient(const FaceArg & face,
230  const StateArg &) const override final;
231  virtual DotType evaluateDot(const ElemArg & elem, const StateArg &) const override final;
232 
237 
248 
249 public:
250  // *********************************************************************************
251  // *********************************************************************************
252  // The following functions are separated here because they are not essential for the
253  // solver but are necessary to interface with the auxiliary and postprocessor
254  // systems.
255  // *********************************************************************************
256  // *********************************************************************************
257 
258  virtual void setDofValue(const DofValue & /*value*/, unsigned int /*index*/) override;
259 
260  virtual void getDofIndices(const Elem * elem,
261  std::vector<dof_id_type> & dof_indices) const override;
262 
263  virtual void setDofValues(const DenseVector<DofValue> & values) override;
264 
265  virtual void clearDofIndices() override;
266 
267  virtual unsigned int numberOfDofs() const override final { return 1; }
268  virtual unsigned int numberOfDofsNeighbor() override final { return 1; }
269 
270  virtual unsigned int oldestSolutionStateRequested() const override final;
271 
272  virtual void clearAllDofIndices() override final;
273 
274  [[noreturn]] virtual const std::vector<dof_id_type> & dofIndicesLower() const override final;
275  [[noreturn]] virtual const FieldVariablePhiValue & phiLower() const override;
276 
277  // Overriding these to make sure nothing happens during residual/jacobian setup.
278  // The only time this can actually happen is when residual setup is called on the auxiliary
279  // system.
280  virtual void residualSetup() override {}
281  virtual void jacobianSetup() override {}
282 
283  virtual libMesh::FEContinuity getContinuity() const override
284  {
285  return _element_data->getContinuity();
286  };
287 
288  virtual void setNodalValue(const OutputType & value, unsigned int idx = 0) override;
289 
290  [[noreturn]] virtual const DofValues & nodalVectorTagValue(TagID) const override;
291 
292  virtual const std::vector<dof_id_type> & dofIndices() const final;
293  virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const final;
294 
295  virtual void prepare() override final {}
296  virtual void prepareNeighbor() override final {}
297  virtual void prepareAux() override final {}
298  virtual void reinitNode() override final {}
299  virtual void reinitNodes(const std::vector<dof_id_type> & /*nodes*/) override final {}
300  virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & /*nodes*/) override final {}
301  virtual void reinitAux() override final {}
302  virtual void reinitAuxNeighbor() override final {}
303  virtual void prepareLowerD() override final {}
304 
305  virtual void computeElemValuesFace() override;
306  virtual void computeNeighborValuesFace() override;
307  virtual void computeNeighborValues() override;
308  virtual void computeLowerDValues() override final;
309 
310  virtual void computeNodalNeighborValues() override final;
311  virtual void computeNodalValues() override final;
312 
313  virtual void computeElemValues() override;
314  virtual void computeFaceValues(const FaceInfo & /*fi*/) override {}
315 
316  virtual void setLowerDofValues(const DenseVector<DofValue> & values) override;
317 
318  virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) override;
319  virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) override;
320  virtual void add(libMesh::NumericVector<libMesh::Number> & vector) override;
321 
322  virtual void setActiveTags(const std::set<TagID> & vtags) override;
323 
324  [[noreturn]] virtual const MooseArray<OutputType> & nodalValueArray() const override;
325  [[noreturn]] virtual const MooseArray<OutputType> & nodalValueOldArray() const override;
326  [[noreturn]] virtual const MooseArray<OutputType> & nodalValueOlderArray() const override;
327 
328  virtual const FieldVariablePhiValue & phi() const override final { return _phi; }
329  virtual const FieldVariablePhiGradient & gradPhi() const override final { return _grad_phi; }
330  [[noreturn]] virtual const FieldVariablePhiSecond & secondPhi() const override final;
331  [[noreturn]] const FieldVariablePhiValue & curlPhi() const override final;
332  [[noreturn]] const FieldVariablePhiDivergence & divPhi() const override final;
333 
334  virtual const FieldVariablePhiValue & phiFace() const override final { return _phi_face; }
335  virtual const FieldVariablePhiGradient & gradPhiFace() const override final
336  {
337  return _grad_phi_face;
338  }
339  [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiFace() const override final;
340 
341  virtual const FieldVariablePhiValue & phiFaceNeighbor() const override final
342  {
343  return _phi_face_neighbor;
344  }
345  virtual const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
346  {
348  }
349  [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final;
350 
351  virtual const FieldVariablePhiValue & phiNeighbor() const override final { return _phi_neighbor; }
352  virtual const FieldVariablePhiGradient & gradPhiNeighbor() const override final
353  {
354  return _grad_phi_neighbor;
355  }
356  [[noreturn]] virtual const FieldVariablePhiSecond & secondPhiNeighbor() const override final;
357 
358  virtual const FieldVariableValue & vectorTagValue(TagID tag) const override;
359  virtual const DofValues & vectorTagDofValue(TagID tag) const override;
360  [[noreturn]] virtual const DofValues & nodalMatrixTagValue(TagID tag) const override;
361  virtual const FieldVariableValue & matrixTagValue(TagID tag) const override;
362 
363  virtual const FieldVariableValue & sln() const override;
364  virtual const FieldVariableValue & slnOld() const override;
365  virtual const FieldVariableValue & slnOlder() const override;
366  virtual const FieldVariableGradient & gradSln() const override;
367  virtual const FieldVariableGradient & gradSlnOld() const override;
368  virtual const FieldVariableValue & slnNeighbor() const override;
369  virtual const FieldVariableValue & slnOldNeighbor() const override;
370  virtual const FieldVariableGradient & gradSlnNeighbor() const override;
371  virtual const FieldVariableGradient & gradSlnOldNeighbor() const override;
372 
373  [[noreturn]] virtual const ADTemplateVariableSecond<OutputType> & adSecondSln() const override;
374  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDot() const override;
375  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDotDot() const override;
376  [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & adGradSlnDot() const override;
377  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adSlnNeighbor() const override;
378  [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> &
379  adGradSlnNeighbor() const override;
380  [[noreturn]] virtual const ADTemplateVariableSecond<OutputType> &
381  adSecondSlnNeighbor() const override;
382  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adUDotNeighbor() const override;
383  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> &
384  adUDotDotNeighbor() const override;
385  [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> &
386  adGradSlnNeighborDot() const override;
387  [[noreturn]] virtual const ADTemplateVariableValue<OutputType> & adSln() const override;
388  [[noreturn]] virtual const ADTemplateVariableGradient<OutputType> & adGradSln() const override;
389  [[noreturn]] virtual const ADTemplateVariableCurl<OutputType> & adCurlSln() const override;
390  [[noreturn]] virtual const ADTemplateVariableCurl<OutputType> &
391  adCurlSlnNeighbor() const override;
392 
393  virtual const DofValues & dofValues() const override;
394  virtual const DofValues & dofValuesOld() const override;
395 
396  virtual const DofValues & dofValuesOlder() const override;
397  virtual const DofValues & dofValuesPreviousNL() const override;
398  virtual const DofValues & dofValuesNeighbor() const override;
399  virtual const DofValues & dofValuesOldNeighbor() const override;
400  virtual const DofValues & dofValuesOlderNeighbor() const override;
401  virtual const DofValues & dofValuesPreviousNLNeighbor() const override;
402  [[noreturn]] virtual const DofValues & dofValuesDot() const override;
403  [[noreturn]] virtual const DofValues & dofValuesDotNeighbor() const override;
404  [[noreturn]] virtual const DofValues & dofValuesDotOld() const override;
405  [[noreturn]] virtual const DofValues & dofValuesDotOldNeighbor() const override;
406  [[noreturn]] virtual const DofValues & dofValuesDotDot() const override;
407  [[noreturn]] virtual const DofValues & dofValuesDotDotNeighbor() const override;
408  [[noreturn]] virtual const DofValues & dofValuesDotDotOld() const override;
409  [[noreturn]] virtual const DofValues & dofValuesDotDotOldNeighbor() const override;
410  [[noreturn]] virtual const MooseArray<libMesh::Number> & dofValuesDuDotDu() const override;
411  [[noreturn]] virtual const MooseArray<libMesh::Number> &
412  dofValuesDuDotDuNeighbor() const override;
413  [[noreturn]] virtual const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const override;
414  [[noreturn]] virtual const MooseArray<libMesh::Number> &
415  dofValuesDuDotDotDuNeighbor() const override;
416 
417  [[noreturn]] virtual const ADDofValues & adDofValues() const override;
418  [[noreturn]] virtual const ADDofValues & adDofValuesNeighbor() const override;
419  [[noreturn]] virtual const ADDofValues & adDofValuesDot() const override;
420  [[noreturn]] virtual const dof_id_type & nodalDofIndex() const override final;
421  [[noreturn]] virtual const dof_id_type & nodalDofIndexNeighbor() const override final;
422 
423  virtual std::size_t phiSize() const override final { return _phi.size(); }
424  virtual std::size_t phiFaceSize() const override final { return _phi_face.size(); }
425  virtual std::size_t phiNeighborSize() const override final { return _phi_neighbor.size(); }
426  virtual std::size_t phiFaceNeighborSize() const override final
427  {
428  return _phi_face_neighbor.size();
429  }
430  [[noreturn]] virtual std::size_t phiLowerSize() const override final;
431 };
432 
433 template <typename OutputType>
434 typename MooseLinearVariableFV<OutputType>::ValueType
435 MooseLinearVariableFV<OutputType>::evaluate(const ElemArg & elem_arg, const StateArg & state) const
436 {
437  const auto & elem_info = this->_mesh.elemInfo(elem_arg.elem->id());
438  return getElemValue(elem_info, state);
439 }
440 
441 template <typename OutputType>
444  const StateArg & state) const
445 {
446  const auto & elem_info = this->_mesh.elemInfo(elem_point.elem->id());
447  return getElemValue(elem_info, state);
448 }
449 
450 template <typename OutputType>
453 {
454  const auto & elem_info = this->_mesh.elemInfo(elem_qp.elem->id());
455  return getElemValue(elem_info, state);
456 }
457 
458 template <typename OutputType>
461  const StateArg & state) const
462 {
463  return (*this)(ElemPointArg{elem_side_qp.elem, elem_side_qp.point, false}, state);
464 }
465 
466 template <typename OutputType>
469  const StateArg & /*state*/) const
470 {
471  const auto & elem_info = this->_mesh.elemInfo(qp_arg.elem->id());
472  return gradSln(elem_info);
473 }
474 
475 template <typename OutputType>
478  const StateArg & /*state*/) const
479 {
480  const auto & elem_info = this->_mesh.elemInfo(elem_arg.elem->id());
481  return gradSln(elem_info);
482 }
483 
484 template <typename OutputType>
487  const StateArg & state) const
488 {
489  mooseAssert(face.fi, "We must have a non-null face information");
490  return gradSln(*face.fi, state);
491 }
492 
493 template <typename OutputType>
494 void
496 {
497  mooseError("MooseLinearVariableFV does not support time integration at the moment! The variable "
498  "which is causing the issue: ",
499  this->name());
500 }
501 
502 template <typename OutputType>
503 void
505 {
506  mooseError("Lower dimensional element support not implemented for finite volume variables!The "
507  "variable which is causing the issue: ",
508  this->name());
509 }
510 
511 template <typename OutputType>
512 void
514 {
515  mooseError("FV variables don't support nodal variable treatment! The variable which is causing "
516  "the issue: ",
517  this->name());
518 }
519 
520 template <typename OutputType>
521 void
523 {
524  mooseError("Linear FV variable does not support automatic differentiation, the variable which is "
525  "attempting it is: ",
526  this->name());
527 }
528 
529 // Declare all the specializations, as the template specialization declarations below must know
530 template <>
531 ADReal MooseLinearVariableFV<Real>::evaluateDot(const ElemArg & elem, const StateArg & state) const;
532 
533 // Prevent implicit instantiation in other translation units where these classes are used
534 extern template class MooseLinearVariableFV<Real>;
virtual const FieldVariablePhiGradient & gradPhiFace() const override final
Return the gradients of the variable&#39;s shape functions on an element face.
virtual void clearAllDofIndices() override final
typename OutputTools< typename Moose::ADType< T >::type >::VariableSecond ADTemplateVariableSecond
Definition: MooseTypes.h:654
virtual const MooseArray< libMesh::Number > & dofValuesDuDotDotDuNeighbor() const override
virtual void setActiveTags(const std::set< TagID > &vtags) override
Set the active vector tags.
typename OutputTools< typename Moose::ADType< T >::type >::VariableCurl ADTemplateVariableCurl
Definition: MooseTypes.h:656
MooseLinearVariableFV< Real > MooseLinearVariableFVReal
virtual bool isExtrapolatedBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const override
Returns whether this (sided) face is an extrapolated boundary face for this functor.
virtual void prepare() override final
Prepare the elemental degrees of freedom.
virtual void computeElemValuesFace() override
Compute values at facial quadrature points.
virtual const DofValues & dofValuesDotDot() const override
Base class for boundary conditions for linear FV systems.
virtual const MooseArray< OutputType > & nodalValueOlderArray() const override
virtual void reinitAuxNeighbor() override final
virtual const ADTemplateVariableCurl< OutputType > & adCurlSln() const override
AD curl solution getter.
virtual DotType evaluateDot(const ElemArg &elem, const StateArg &) const override final
Evaluate the functor time derivative with a given element.
void adError() const
Throw an error when somebody wants to use this variable with automatic differentiation.
std::unordered_map< BoundaryID, LinearFVBoundaryCondition * > _boundary_id_to_bc
Map for easily accessing the boundary conditions based on the boundary IDs.
virtual const FieldVariableValue & sln() const override
virtual bool isDirichletBoundaryFace(const FaceInfo &fi) const
If the variable has a dirichlet boundary condition at face described by fi .
virtual const DofValues & dofValuesDotOld() const override
unsigned int TagID
Definition: MooseTypes.h:238
const FieldVariablePhiGradient & _grad_phi_face
Class for stuff related to variables.
const FieldVariablePhiGradient & _grad_phi_face_neighbor
virtual const MooseArray< OutputType > & nodalValueArray() const override
Methods for retrieving values of variables at the nodes in a MooseArray for AuxKernelBase.
const FieldVariablePhiGradient & _grad_phi_neighbor
virtual const Elem *const & currentElem() const override
Current element this variable is evaluated at.
virtual const DofValues & dofValuesDotDotNeighbor() const override
virtual const ADTemplateVariableGradient< OutputType > & adGradSlnDot() const override
AD grad of time derivative solution getter.
virtual 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...
virtual bool isNodal() const override final
Is this variable nodal.
typename MooseVariableField< OutputType >::FieldVariablePhiSecond FieldVariablePhiSecond
const FieldVariablePhiDivergence & divPhi() const override final
Divergence of the shape functions.
virtual const FieldVariablePhiSecond & secondPhiFaceNeighbor() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s shape functions on a neighboring ele...
virtual const DofValues & dofValues() const override
dof values getters
virtual const ADTemplateVariableValue< OutputType > & adSln() const override
AD solution getter.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
virtual const FieldVariableGradient & gradSlnOldNeighbor() const override
virtual void residualSetup() override
virtual const DofValues & nodalVectorTagValue(TagID) const override
Base class for finite volume Dirichlet boundaray conditions.
const unsigned int _sys_num
Cache the number of the system this variable belongs to.
virtual std::size_t phiFaceSize() const override final
Return phiFace size.
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:648
virtual const DofValues & dofValuesPreviousNL() const override
virtual const DofValues & nodalMatrixTagValue(TagID tag) const override
virtual const DofValues & dofValuesPreviousNLNeighbor() const override
virtual const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void setLowerDofValues(const DenseVector< DofValue > &values) override
virtual void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const override
virtual const DofValues & dofValuesDotDotOldNeighbor() const override
virtual const FieldVariablePhiValue & phi() const override final
Return the variable&#39;s elemental shape functions.
virtual const FieldVariablePhiGradient & gradPhi() const override final
Return the gradients of the variable&#39;s elemental shape functions.
RealVectorValue _cell_gradient
Temporary storage for the cell gradient to avoid unnecessary allocations.
virtual void prepareNeighbor() override final
Prepare the neighbor element degrees of freedom.
virtual void setNodalValue(const OutputType &value, unsigned int idx=0) override
LinearFVBoundaryCondition * getBoundaryCondition(const BoundaryID bd_id) const
Get the boundary condition object which corresponds to the given boundary ID.
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
virtual const FieldVariablePhiSecond & secondPhi() const override final
Return the rank-2 tensor of second derivatives of the variable&#39;s elemental shape functions.
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...
virtual void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
const FieldVariablePhiValue & _phi_face_neighbor
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:42
virtual const DofValues & dofValuesDotDotOld() const override
virtual bool supportsElemSideQpArg() const override final
Whether this functor supports evaluation with ElemSideQpArg.
virtual void clearDofIndices() override
Clear out the dof indices.
virtual const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
virtual void computeElemValues() override
Compute values at interior quadrature points.
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
virtual const ADTemplateVariableCurl< OutputType > & adCurlSlnNeighbor() const override
AD curl neighbor solution getter.
void lowerDError() const
Throw and error when somebody requests lower-dimensional data from this variable. ...
virtual const MooseArray< libMesh::Number > & dofValuesDuDotDuNeighbor() const override
virtual 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...
virtual const ADTemplateVariableValue< OutputType > & adUDotDot() const override
AD second time derivative getter.
virtual void prepareLowerD() override final
Prepare a lower dimensional element&#39;s degrees of freedom.
libMesh::TensorTools::DecrementRank< OutputShape >::type OutputShapeDivergence
virtual const ADTemplateVariableGradient< OutputType > & adGradSlnNeighbor() const override
AD grad neighbor solution getter.
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
virtual const ADDofValues & adDofValuesNeighbor() const override
Return the AD neighbor dof values.
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
virtual bool supportsFaceArg() const override final
Whether this functor supports evaluation with FaceArg.
virtual bool isFV() const override
virtual const DofValues & dofValuesDotNeighbor() const override
virtual const FieldVariableValue & slnOld() const override
typename MooseVariableField< OutputType >::FieldVariablePhiValue FieldVariablePhiValue
void timeIntegratorError() const
Throw an error when somebody requests time-related data from this variable.
dof_id_type id() const
virtual std::size_t phiSize() const override final
Return phi size.
A structure defining a "face" evaluation calling argument for Moose functors.
const std::vector< std::unique_ptr< libMesh::NumericVector< libMesh::Number > > > & _grad_container
Pointer to the cell gradients which are stored on the linear system.
typename MooseVariableField< OutputType >::FieldVariablePhiGradient FieldVariablePhiGradient
virtual const DofValues & vectorTagDofValue(TagID tag) const override
virtual void setDofValue(const DofValue &, unsigned int) override
virtual void insertLower(libMesh::NumericVector< libMesh::Number > &vector) override
Insert the currently cached degree of freedom values for a lower-dimensional element into the provide...
const FaceInfo * fi
a face information object which defines our location in space
virtual const FieldVariablePhiValue & phiLower() const override
Return the variable&#39;s shape functions on a lower-dimensional element.
boundary_id_type BoundaryID
Point point
The physical location of the quadrature point.
typename MooseVariableField< OutputType >::DofValues DofValues
virtual void reinitNode() override final
const FieldVariablePhiGradient & _grad_phi
virtual const FieldVariableValue & matrixTagValue(TagID tag) const override
const libMesh::Elem * elem
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:651
const libMesh::Elem * elem
virtual void reinitNodes(const std::vector< dof_id_type > &) override final
Provides an interface for computing residual contributions from finite volume numerical fluxes comput...
Definition: FVFluxBC.h:23
virtual const DofValues & dofValuesNeighbor() const override
std::unique_ptr< MooseVariableDataLinearFV< OutputType > > _element_data
Holder for all the data associated with the "main" element.
virtual const ADTemplateVariableGradient< OutputType > & adGradSlnNeighborDot() const override
AD grad of time derivative neighbor solution getter.
virtual const DofValues & dofValuesOlder() const override
virtual unsigned int numberOfDofsNeighbor() override final
A structure that is used to evaluate Moose functors logically at an element/cell center.
typename MooseVariableField< OutputType >::FieldVariableGradient FieldVariableGradient
MooseLinearVariableFV(const InputParameters &parameters)
Argument for requesting functor evaluation at a quadrature point location in an element.
virtual const dof_id_type & nodalDofIndex() const override final
virtual const MooseArray< libMesh::Number > & dofValuesDuDotDu() const override
const libMesh::Elem * elem
The element.
virtual const ADTemplateVariableSecond< OutputType > & adSecondSln() const override
AD second solution getter.
const FieldVariablePhiValue & curlPhi() const override final
Curl of the shape functions.
void nodalError() const
Throw an error when somebody wants to use this variable as a nodal variable.
virtual void jacobianSetup() override
virtual void prepareAux() override final
virtual void reinitNodesNeighbor(const std::vector< dof_id_type > &) override final
virtual bool hasDoFsOnNodes() const override final
Does this variable have DoFs on nodes.
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 bool computingSecond() const override final
Whether or not this variable is computing any second derivatives.
virtual const FieldVariablePhiGradient & gradPhiFaceNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element face.
virtual const FieldVariableValue & slnNeighbor() const override
virtual bool isNodalDefined() const override final
Is this variable defined at nodes.
virtual const ADDofValues & adDofValuesDot() const override
Return the AD time derivatives at dofs.
virtual const FieldVariableGradient & gradSlnNeighbor() const override
neighbor solution gradients
virtual const FieldVariableGradient & gradSlnOld() const override
virtual const FieldVariablePhiValue & phiNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element.
virtual const MooseArray< OutputType > & nodalValueOldArray() const override
MooseMesh & _mesh
mesh the variable is active in
typename MooseVariableDataBase< Real >::ADDofValue ADDofValue
virtual const DofValues & dofValuesDotOldNeighbor() const override
virtual void computeNodalNeighborValues() override final
Compute nodal values of this variable in the neighbor.
virtual const DofValues & dofValuesOldNeighbor() const override
virtual void reinitAux() override final
virtual void add(libMesh::NumericVector< libMesh::Number > &vector) override
Add the currently cached degree of freedom values into the provided vector.
typename MooseVariableField< OutputType >::ADDofValues ADDofValues
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeNeighborValuesFace() override
Compute values at facial quadrature points for the neighbor.
virtual bool usesSecondPhiNeighbor() const override final
Whether or not this variable is actually using the shape function second derivatives.
virtual const FieldVariableValue & vectorTagValue(TagID tag) const override
tag values getters
const FieldVariablePhiValue & _phi_face
virtual void sizeMatrixTagData() override
Size data structures related to matrix tagging.
std::unique_ptr< MooseVariableDataLinearFV< OutputType > > _neighbor_data
Holder for all the data associated with the "neighbor" element.
virtual void computeNodalValues() override final
Compute nodal values of this variable.
typename MooseVariableDataBase< Real >::ADDofValues ADDofValues
const std::unordered_map< BoundaryID, LinearFVBoundaryCondition * > & getBoundaryConditionMap()
static InputParameters validParams()
virtual std::size_t phiFaceNeighborSize() const override final
Return phiFaceNeighbor size.
virtual const dof_id_type & nodalDofIndexNeighbor() const override final
typename MooseVariableField< OutputType >::FieldVariablePhiDivergence FieldVariablePhiDivergence
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
Get the solution value for the provided element and seed the derivative for the corresponding dof ind...
virtual void computeLowerDValues() override final
compute values at quadrature points on the lower dimensional element
void cacheBoundaryBCMap()
Setup the boundary to Dirichlet BC map.
virtual std::size_t phiLowerSize() const override final
Return the number of shape functions on the lower dimensional element for this variable.
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
const FieldVariablePhiValue & _phi
Shape functions, only used when we are postprocessing or using this variable in an auxiliary system...
typename MooseVariableField< OutputType >::FieldVariableValue FieldVariableValue
virtual void insert(libMesh::NumericVector< libMesh::Number > &vector) override
Insert the currently cached degree of freedom values into the provided vector.
virtual const ADTemplateVariableValue< OutputType > & adUDotNeighbor() const override
AD neighbor time derivative getter.
State argument for evaluating functors.
virtual std::size_t phiNeighborSize() const override final
Return phiNeighbor size.
typename MooseVariableDataBase< Real >::DofValues DofValues
virtual void computeFaceValues(const FaceInfo &) override
Compute values at face quadrature points for the element+neighbor (both sides of the face)...
virtual const ADTemplateVariableValue< OutputType > & adSlnNeighbor() const override
AD neighbor solution getter.
virtual void prepareIC() override
Prepare the initial condition.
virtual bool needsGradientVectorStorage() const override
Check if cell gradient computations were requested for this variable.
virtual bool computingCurl() const override final
Whether or not this variable is computing any curl quantities.
libMesh::TensorTools::DecrementRank< Real >::type OutputDivergence
virtual const ADTemplateVariableSecond< OutputType > & adSecondSlnNeighbor() const override
AD second neighbor solution getter.
virtual void computeNeighborValues() override
Compute values at quadrature points for the neighbor.
Class used for caching additional information for elements such as the volume and centroid...
Definition: ElemInfo.h:25
virtual const DofValues & dofValuesOlderNeighbor() const override
virtual libMesh::FEContinuity getContinuity() const override
Return the continuity of this variable.
const libMesh::NumericVector< libMesh::Number > *const & _solution
The current (ghosted) solution.
Real Number
virtual const FieldVariablePhiGradient & gradPhiNeighbor() const override final
Return the gradients of the variable&#39;s shape functions on a neighboring element.
void computeCellGradients()
Switch to request cell gradient computations.
virtual const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const override
virtual const ADDofValues & adDofValues() const override
Return the AD dof values.
virtual unsigned int oldestSolutionStateRequested() const override final
The oldest solution state that is requested for this variable (0 = current, 1 = old, 2 = older, etc).
typename MooseVariableDataBase< Real >::DofValue DofValue
virtual void setDofValues(const DenseVector< DofValue > &values) override
const libMesh::Elem * elem
The element.
virtual const DofValues & dofValuesOld() const override
const FieldVariablePhiValue & _phi_neighbor
This class provides variable solution interface for linear finite volume problems.
Definition: FVUtils.h:24
virtual unsigned int numberOfDofs() const override final
Get the number of local DoFs.
virtual const FieldVariableValue & slnOlder() const override
virtual const ADTemplateVariableValue< OutputType > & adUDot() const override
AD time derivative getter.
virtual bool computingDiv() const override final
Whether or not this variable is computing any divergence quantities.
virtual const FieldVariableGradient & gradSln() const override
element gradients
virtual const ADTemplateVariableGradient< OutputType > & adGradSln() const override
AD grad solution getter.
virtual GradientType evaluateGradient(const ElemQpArg &qp_arg, const StateArg &) const override final
Argument for requesting functor evaluation at quadrature point locations on an element side...
virtual const FieldVariablePhiValue & phiFaceNeighbor() const override final
Return the variable&#39;s shape functions on a neighboring element face.
const ElemInfo & elemInfo(const dof_id_type id) const
Accessor for the elemInfo object for a given element ID.
Definition: MooseMesh.C:3951
Moose::ShapeType< Real >::type OutputShape
virtual const DofValues & dofValuesDot() const override
uint8_t dof_id_type
bool _needs_cell_gradients
Boolean to check if this variable needs gradient computations.
virtual const FieldVariableValue & slnOldNeighbor() const override
virtual ValueType evaluate(const ElemArg &elem, const StateArg &) const override final
Evaluate the functor with a given element.
virtual const ADTemplateVariableValue< OutputType > & adUDotDotNeighbor() const override
AD neighbor second time derivative getter.