https://mooseframework.inl.gov
MooseVariableFieldBase.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 "MooseVariableBase.h"
13 
14 namespace libMesh
15 {
16 template <typename>
17 class DenseVector;
18 template <typename>
19 class NumericVector;
20 class Point;
21 }
22 
23 class FaceInfo;
24 
30 {
31 public:
33 
35 
40  virtual void clearDofIndices() = 0;
41 
45  virtual void prepare() = 0;
46 
50  virtual void prepareNeighbor() = 0;
51 
55  virtual void prepareLowerD() = 0;
56 
57  virtual void prepareAux() = 0;
58 
59  virtual void reinitNode() = 0;
60  virtual void reinitAux() = 0;
61  virtual void reinitAuxNeighbor() = 0;
62 
63  virtual void reinitNodes(const std::vector<dof_id_type> & nodes) = 0;
64  virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & nodes) = 0;
65 
69  virtual Moose::VarFieldType fieldType() const = 0;
70 
74  const std::string & componentName(const unsigned int comp) const;
75 
79  virtual bool isVector() const = 0;
80 
81  virtual bool isFV() const { return false; }
82 
87  virtual bool isNodalDefined() const = 0;
88 
89  virtual const dof_id_type & nodalDofIndex() const = 0;
90  virtual const dof_id_type & nodalDofIndexNeighbor() const = 0;
91 
95  virtual const Elem * const & currentElem() const = 0;
96 
100  const std::set<SubdomainID> & activeSubdomains() const;
101 
107  bool activeOnSubdomain(SubdomainID subdomain) const;
108 
114  bool activeOnSubdomains(const std::set<SubdomainID> & subdomains) const;
115 
120  virtual bool needsGradientVectorStorage() const { return false; }
121 
125  virtual void prepareIC() = 0;
126 
131  virtual void computeFaceValues(const FaceInfo & /*fi*/) { mooseError("not implemented"); }
132 
136  virtual void computeElemValues() = 0;
140  virtual void computeElemValuesFace() = 0;
144  virtual void computeNeighborValuesFace() = 0;
148  virtual void computeNeighborValues() = 0;
152  virtual void computeLowerDValues() = 0;
156  virtual void computeNodalNeighborValues() = 0;
160  virtual void computeNodalValues() = 0;
161 
166  virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const = 0;
167 
173  virtual const std::vector<dof_id_type> & dofIndicesLower() const = 0;
174 
175  virtual unsigned int numberOfDofsNeighbor() = 0;
176 
180  virtual void insert(libMesh::NumericVector<libMesh::Number> & vector) = 0;
181 
186  virtual void insertLower(libMesh::NumericVector<libMesh::Number> & vector) = 0;
187 
191  virtual void add(libMesh::NumericVector<libMesh::Number> & vector) = 0;
192 
196  virtual std::size_t phiSize() const = 0;
200  virtual std::size_t phiFaceSize() const = 0;
204  virtual std::size_t phiNeighborSize() const = 0;
208  virtual std::size_t phiFaceNeighborSize() const = 0;
212  virtual std::size_t phiLowerSize() const = 0;
213 
218  virtual unsigned int oldestSolutionStateRequested() const = 0;
219 };
220 
221 #define usingMooseVariableFieldBaseMembers usingMooseVariableBaseMembers
virtual unsigned int numberOfDofsNeighbor()=0
VarFieldType
Definition: MooseTypes.h:722
virtual void clearDofIndices()=0
Clear out the dof indices.
static InputParameters validParams()
virtual bool needsGradientVectorStorage() const
Check if this variable needs a raw vector of gradients at dof-values.
virtual void prepareNeighbor()=0
Prepare the neighbor element degrees of freedom.
virtual void reinitAuxNeighbor()=0
virtual void computeNodalValues()=0
Compute nodal values of this variable.
virtual void computeNeighborValuesFace()=0
Compute values at facial quadrature points for the neighbor.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:127
virtual std::size_t phiLowerSize() const =0
Return the number of shape functions on the lower dimensional element for this variable.
virtual bool isFV() const
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void prepareIC()=0
Prepare the initial condition.
virtual void reinitAux()=0
This class provides an interface for common operations on field variables of both FE and FV types wit...
virtual const dof_id_type & nodalDofIndex() const =0
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
virtual void computeFaceValues(const FaceInfo &)
Compute values at face quadrature points for the element+neighbor (both sides of the face)...
virtual bool isVector() const =0
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
virtual void reinitNodesNeighbor(const std::vector< dof_id_type > &nodes)=0
virtual std::size_t phiFaceNeighborSize() const =0
Return phiFaceNeighbor size.
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const =0
Get neighbor DOF indices for currently selected element.
bool activeOnSubdomains(const std::set< SubdomainID > &subdomains) const
Is the variable active on the subdomains?
virtual std::size_t phiFaceSize() const =0
Return phiFace size.
virtual void computeLowerDValues()=0
compute values at quadrature points on the lower dimensional element
virtual void add(libMesh::NumericVector< libMesh::Number > &vector)=0
Add the currently cached degree of freedom values into the provided vector.
virtual void reinitNodes(const std::vector< dof_id_type > &nodes)=0
virtual void reinitNode()=0
virtual std::size_t phiNeighborSize() const =0
Return phiNeighbor size.
virtual const Elem *const & currentElem() const =0
Current element this variable is evaluated at.
virtual void prepare()=0
Prepare the elemental degrees of freedom.
virtual unsigned int oldestSolutionStateRequested() const =0
The oldest solution state that is requested for this variable (0 = current, 1 = old, 2 = older, etc).
virtual void insert(libMesh::NumericVector< libMesh::Number > &vector)=0
Insert the currently cached degree of freedom values into the provided vector.
virtual void insertLower(libMesh::NumericVector< libMesh::Number > &vector)=0
Insert the currently cached degree of freedom values for a lower-dimensional element into the provide...
MooseVariableFieldBase(const InputParameters &parameters)
const std::string & componentName(const unsigned int comp) const
Get the variable name of a component in libMesh.
virtual bool isNodalDefined() const =0
Is this variable defined at nodes.
virtual void computeElemValues()=0
Compute values at interior quadrature points.
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
virtual void computeElemValuesFace()=0
Compute values at facial quadrature points.
virtual void prepareAux()=0
virtual void computeNeighborValues()=0
Compute values at quadrature points for the neighbor.
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:267
virtual Moose::VarFieldType fieldType() const =0
Field type of this variable.
virtual std::size_t phiSize() const =0
Return phi size.
bool activeOnSubdomain(SubdomainID subdomain) const
Is the variable active on the subdomain?
virtual void computeNodalNeighborValues()=0
Compute nodal values of this variable in the neighbor.
virtual const std::vector< dof_id_type > & dofIndicesLower() const =0
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
virtual const dof_id_type & nodalDofIndexNeighbor() const =0
virtual void prepareLowerD()=0
Prepare a lower dimensional element&#39;s degrees of freedom.
Base variable class.
uint8_t dof_id_type