www.mooseframework.org
MooseVariableFEBase.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 Assembly;
24 
26 {
27 public:
28  MooseVariableFEBase(unsigned int var_num,
29  const FEType & fe_type,
30  SystemBase & sys,
31  Moose::VarKindType var_kind,
32  THREAD_ID tid);
33 
38  virtual void clearDofIndices() = 0;
39 
43  virtual void prepare() = 0;
44 
48  virtual void prepareNeighbor() = 0;
49 
53  virtual void prepareLowerD() = 0;
54 
55  virtual void prepareAux() = 0;
56 
57  virtual void reinitNode() = 0;
58  virtual void reinitAux() = 0;
59  virtual void reinitAuxNeighbor() = 0;
60 
61  virtual void reinitNodes(const std::vector<dof_id_type> & nodes) = 0;
62  virtual void reinitNodesNeighbor(const std::vector<dof_id_type> & nodes) = 0;
63 
68  virtual bool isNodal() const = 0;
69 
73  virtual bool isVector() const = 0;
74 
75  virtual const dof_id_type & nodalDofIndex() const = 0;
76  virtual const dof_id_type & nodalDofIndexNeighbor() const = 0;
77 
81  virtual const Elem * const & currentElem() const = 0;
82 
86  virtual const std::set<SubdomainID> & activeSubdomains() const = 0;
92  virtual bool activeOnSubdomain(SubdomainID subdomain) const = 0;
93 
97  virtual void prepareIC() = 0;
98 
102  virtual void computeElemValues() = 0;
106  virtual void computeElemValuesFace() = 0;
110  virtual void computeNeighborValuesFace() = 0;
114  virtual void computeNeighborValues() = 0;
118  virtual void computeLowerDValues() = 0;
122  virtual void computeNodalNeighborValues() = 0;
126  virtual void computeNodalValues() = 0;
130  virtual void setDofValues(const DenseVector<Number> & value) = 0;
134  virtual Number getNodalValue(const Node & node) = 0;
138  virtual Number getNodalValueOld(const Node & node) = 0;
142  virtual Number getNodalValueOlder(const Node & node) = 0;
149  virtual Number getElementalValue(const Elem * elem, unsigned int idx = 0) const = 0;
156  virtual Number getElementalValueOld(const Elem * elem, unsigned int idx = 0) const = 0;
163  virtual Number getElementalValueOlder(const Elem * elem, unsigned int idx = 0) const = 0;
164 
165  virtual void getDofIndices(const Elem * elem, std::vector<dof_id_type> & dof_indices) const = 0;
170  virtual const std::vector<dof_id_type> & dofIndicesNeighbor() const = 0;
171 
177  virtual const std::vector<dof_id_type> & dofIndicesLower() const = 0;
178 
179  virtual unsigned int numberOfDofsNeighbor() = 0;
180 
181  virtual void insert(NumericVector<Number> & residual) = 0;
182  virtual void add(NumericVector<Number> & residual) = 0;
183 
187  virtual const MooseArray<Number> & dofValue() = 0;
191  virtual const MooseArray<Number> & dofValues() = 0;
195  virtual const MooseArray<Number> & dofValuesOld() = 0;
199  virtual const MooseArray<Number> & dofValuesOlder() = 0;
203  virtual const MooseArray<Number> & dofValuesPreviousNL() = 0;
207  virtual const MooseArray<Number> & dofValuesNeighbor() = 0;
211  virtual const MooseArray<Number> & dofValuesOldNeighbor() = 0;
215  virtual const MooseArray<Number> & dofValuesOlderNeighbor() = 0;
219  virtual const MooseArray<Number> & dofValuesPreviousNLNeighbor() = 0;
223  virtual const MooseArray<Number> & dofValuesDot() = 0;
227  virtual const MooseArray<Number> & dofValuesDotNeighbor() = 0;
231  virtual const MooseArray<Number> & dofValuesDotDot() = 0;
235  virtual const MooseArray<Number> & dofValuesDotDotNeighbor() = 0;
239  virtual const MooseArray<Number> & dofValuesDotOld() = 0;
243  virtual const MooseArray<Number> & dofValuesDotOldNeighbor() = 0;
247  virtual const MooseArray<Number> & dofValuesDotDotOld() = 0;
251  virtual const MooseArray<Number> & dofValuesDotDotOldNeighbor() = 0;
255  virtual const MooseArray<Number> & dofValuesDuDotDu() = 0;
259  virtual const MooseArray<Number> & dofValuesDuDotDuNeighbor() = 0;
263  virtual const MooseArray<Number> & dofValuesDuDotDotDu() = 0;
267  virtual const MooseArray<Number> & dofValuesDuDotDotDuNeighbor() = 0;
268 
272  virtual size_t phiSize() const = 0;
276  virtual size_t phiFaceSize() const = 0;
280  virtual size_t phiNeighborSize() const = 0;
284  virtual size_t phiFaceNeighborSize() const = 0;
288  virtual size_t phiLowerSize() const = 0;
289 };
290 
virtual unsigned int numberOfDofsNeighbor()=0
virtual Number getElementalValueOlder(const Elem *elem, unsigned int idx=0) const =0
Get the older value of this variable on an element.
virtual void add(NumericVector< Number > &residual)=0
virtual const MooseArray< Number > & dofValuesDotOldNeighbor()=0
Returns old time derivative of neighboring degrees of freedom.
virtual void prepareNeighbor()=0
Prepare the neighbor element degrees of freedom.
virtual void reinitAux()=0
virtual void prepareAux()=0
virtual const MooseArray< Number > & dofValuesDuDotDotDu()=0
Returns derivative of second time derivative of degrees of freedom.
Keeps track of stuff related to assembling.
Definition: Assembly.h:62
virtual void prepare()=0
Prepare the elemental degrees of freedom.
virtual Number getNodalValue(const Node &node)=0
Get the value of this variable at given node.
virtual const MooseArray< Number > & dofValuesDuDotDotDuNeighbor()=0
Returns derivative of second time derivative of neighboring degrees of freedom.
virtual Number getElementalValue(const Elem *elem, unsigned int idx=0) const =0
Get the current value of this variable on an element.
virtual bool activeOnSubdomain(SubdomainID subdomain) const =0
Is the variable active on the subdomain?
virtual const MooseArray< Number > & dofValuesPreviousNLNeighbor()=0
Returns previous nl solution on neighbor element.
virtual const MooseArray< Number > & dofValuesPreviousNL()=0
Returns previous nl solution on element.
virtual void reinitNodesNeighbor(const std::vector< dof_id_type > &nodes)=0
virtual size_t phiLowerSize() const =0
Return the number of shape functions on the lower dimensional element for this variable.
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
virtual const MooseArray< Number > & dofValuesDot()=0
Returns time derivative of degrees of freedom.
Base class for a system (of equations)
Definition: SystemBase.h:92
virtual void prepareIC()=0
Prepare the initial condition.
virtual size_t phiNeighborSize() const =0
Return phiNeighbor size.
virtual const MooseArray< Number > & dofValuesDotNeighbor()=0
Returns time derivative of neighboring degrees of freedom.
virtual const MooseArray< Number > & dofValuesOlderNeighbor()=0
Returns older dof solution on neighbor element.
virtual void computeElemValues()=0
Compute values at interior quadrature points.
virtual size_t phiSize() const =0
Return phi size.
virtual void computeNeighborValues()=0
Compute values at quadrature points for the neighbor.
virtual const dof_id_type & nodalDofIndex() const =0
virtual const std::set< SubdomainID > & activeSubdomains() const =0
The subdomains the variable is active on.
virtual bool isNodal() const =0
Is this variable nodal.
virtual void insert(NumericVector< Number > &residual)=0
virtual void computeNodalNeighborValues()=0
Compute nodal values of this variable in the neighbor.
virtual void computeNeighborValuesFace()=0
Compute values at facial quadrature points for the neighbor.
virtual void computeLowerDValues()=0
compute values at quadrature points on the lower dimensional element
virtual const MooseArray< Number > & dofValues()=0
Returns dof solution on element.
virtual const MooseArray< Number > & dofValuesOld()=0
Returns old dof solution on element.
MooseVariableFEBase(unsigned int var_num, const FEType &fe_type, SystemBase &sys, Moose::VarKindType var_kind, THREAD_ID tid)
virtual const MooseArray< Number > & dofValuesDotDotNeighbor()=0
Returns second time derivative of neighboring degrees of freedom.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:481
virtual void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const =0
virtual void clearDofIndices()=0
Clear out the dof indices.
virtual void reinitNodes(const std::vector< dof_id_type > &nodes)=0
virtual void reinitNode()=0
subdomain_id_type SubdomainID
virtual void computeElemValuesFace()=0
Compute values at facial quadrature points.
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 size_t phiFaceSize() const =0
Return phiFace size.
virtual const MooseArray< Number > & dofValuesDotDotOld()=0
Returns old second time derivative of degrees of freedom.
virtual void prepareLowerD()=0
Prepare a lower dimensional element&#39;s degrees of freedom.
virtual const MooseArray< Number > & dofValuesDotDotOldNeighbor()=0
Returns old second time derivative of neighboring degrees of freedom.
forward declarations
Definition: MooseArray.h:16
virtual size_t phiFaceNeighborSize() const =0
Return phiFaceNeighbor size.
virtual bool isVector() const =0
virtual void reinitAuxNeighbor()=0
virtual Number getElementalValueOld(const Elem *elem, unsigned int idx=0) const =0
Get the old value of this variable on an element.
virtual Number getNodalValueOld(const Node &node)=0
Get the old value of this variable at given node.
virtual const MooseArray< Number > & dofValuesDuDotDu()=0
Returns derivative of time derivative of degrees of freedom.
virtual const Elem *const & currentElem() const =0
Current element this variable is evaluated at.
virtual void computeNodalValues()=0
Compute nodal values of this variable.
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const =0
Get neighbor DOF indices for currently selected element.
virtual const MooseArray< Number > & dofValue()=0
Deprecated method.
virtual const dof_id_type & nodalDofIndexNeighbor() const =0
virtual const MooseArray< Number > & dofValuesOldNeighbor()=0
Returns old dof solution on neighbor element.
SystemBase & sys()
Get the system this variable is part of.
virtual const MooseArray< Number > & dofValuesOlder()=0
Returns older dof solution on element.
virtual const MooseArray< Number > & dofValuesDotDot()=0
Returns second time derivative of degrees of freedom.
virtual const MooseArray< Number > & dofValuesDotOld()=0
Returns old time derivative of degrees of freedom.
virtual void setDofValues(const DenseVector< Number > &value)=0
Set values for this variable to keep everything up to date.
unsigned int THREAD_ID
Definition: MooseTypes.h:161
virtual const MooseArray< Number > & dofValuesDuDotDuNeighbor()=0
Returns derivative of time derivative of neighboring degrees of freedom.
virtual const MooseArray< Number > & dofValuesNeighbor()=0
Returns dof solution on neighbor element.
virtual Number getNodalValueOlder(const Node &node)=0
Get the t-2 value of this variable at given node.