www.mooseframework.org
MooseVariableFE.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 "metaphysicl/numberarray.h"
13 #include "metaphysicl/dualnumber.h"
14 
15 #include "MooseTypes.h"
16 #include "MooseVariableFEBase.h"
17 #include "SubProblem.h"
18 #include "SystemBase.h"
19 #include "MooseMesh.h"
20 #include "MooseVariableData.h"
21 
22 #include "libmesh/numeric_vector.h"
23 #include "libmesh/dof_map.h"
24 #include "libmesh/elem.h"
25 #include "libmesh/quadrature.h"
26 #include "libmesh/dense_vector.h"
27 #include "libmesh/dense_vector.h"
28 
29 class TimeIntegrator;
30 
36 template <typename OutputType>
38 {
39 public:
40  typedef OutputType OutputShape;
41  typedef OutputType OutputValue;
45 
51 
57 
63 
64  MooseVariableFE(unsigned int var_num,
65  const FEType & fe_type,
66  SystemBase & sys,
67  const Assembly & assembly,
68  Moose::VarKindType var_kind,
69  THREAD_ID tid);
70 
71  void clearDofIndices() override;
72 
73  void prepare() override;
74  void prepareNeighbor() override;
75  void prepareLowerD() override;
76 
77  void prepareAux() override;
78 
79  void reinitNode() override;
80  void reinitAux() override;
81  void reinitAuxNeighbor() override;
82 
83  void reinitNodes(const std::vector<dof_id_type> & nodes) override;
84  void reinitNodesNeighbor(const std::vector<dof_id_type> & nodes) override;
85 
91  bool usesPhi() const { return true; }
97  bool usesGradPhi() const { return true; }
103  bool usesPhiNeighbor() const { return true; }
109  bool usesGradPhiNeighbor() const { return true; }
113  bool usesSecondPhi() const;
114 
119  bool usesSecondPhiNeighbor() const;
120 
124  bool computingSecond() const { return usesSecondPhi(); }
125 
129  bool computingCurl() const;
130 
131  const std::set<SubdomainID> & activeSubdomains() const override;
132  bool activeOnSubdomain(SubdomainID subdomain) const override;
133 
134  bool isNodal() const override { return _element_data->isNodal(); }
135  bool isVector() const override;
136  const Node * const & node() const { return _element_data->node(); }
137  const dof_id_type & nodalDofIndex() const override { return _element_data->nodalDofIndex(); }
138  bool isNodalDefined() const;
139 
140  const Node * const & nodeNeighbor() const { return _neighbor_data->node(); }
141  const dof_id_type & nodalDofIndexNeighbor() const override
142  {
143  return _neighbor_data->nodalDofIndex();
144  }
145  bool isNodalNeighborDefined() const;
146 
147  const Elem * const & currentElem() const override { return _element_data->currentElem(); }
148 
152  const unsigned int & currentSide() const { return _element_data->currentSide(); }
153 
157  const Elem * const & neighbor() const { return _neighbor_data->currentElem(); }
158 
159  virtual void prepareIC() override;
160 
161  const FieldVariablePhiValue & phi() const { return _element_data->phi(); }
162  const FieldVariablePhiGradient & gradPhi() const { return _element_data->gradPhi(); }
163  const FieldVariablePhiSecond & secondPhi() const;
164  const FieldVariablePhiCurl & curlPhi() const;
165 
166  const FieldVariablePhiValue & phiFace() const { return _element_data->phiFace(); }
167  const FieldVariablePhiGradient & gradPhiFace() const { return _element_data->gradPhiFace(); }
168  const FieldVariablePhiSecond & secondPhiFace() const;
169  const FieldVariablePhiCurl & curlPhiFace() const;
170 
171  const FieldVariablePhiValue & phiNeighbor() const { return _neighbor_data->phi(); }
172  const FieldVariablePhiGradient & gradPhiNeighbor() const { return _neighbor_data->gradPhi(); }
174  const FieldVariablePhiCurl & curlPhiNeighbor() const;
175 
176  const FieldVariablePhiValue & phiFaceNeighbor() const { return _neighbor_data->phiFace(); }
178  {
179  return _neighbor_data->gradPhiFace();
180  }
183 
184  const FieldVariablePhiValue & phiLower() const { return _lower_data->phi(); }
185  const FieldVariablePhiGradient & gradPhiLower() const { return _lower_data->gradPhi(); }
186 
187  template <ComputeStage compute_stage>
189  {
190  return _element_data->template adGradPhi<compute_stage>();
191  }
192 
193  template <ComputeStage compute_stage>
195  {
196  return _element_data->template adGradPhiFace<compute_stage>();
197  }
198 
199  // damping
200  const FieldVariableValue & increment() const { return _element_data->increment(); }
201 
203  {
204  return _element_data->vectorTagValue(tag);
205  }
207  {
208  return _element_data->matrixTagValue(tag);
209  }
210 
212  const FieldVariableValue & sln() const { return _element_data->sln(Moose::Current); }
213  const FieldVariableValue & slnOld() const { return _element_data->sln(Moose::Old); }
214  const FieldVariableValue & slnOlder() const { return _element_data->sln(Moose::Older); }
216 
218  const FieldVariableGradient & gradSln() const { return _element_data->gradSln(Moose::Current); }
219  const FieldVariableGradient & gradSlnOld() const { return _element_data->gradSln(Moose::Old); }
221  {
222  return _element_data->gradSln(Moose::Older);
223  }
225  {
226  return _element_data->gradSln(Moose::PreviousNL);
227  }
228 
230  const FieldVariableGradient & gradSlnDot() const { return _element_data->gradSlnDot(); }
231  const FieldVariableGradient & gradSlnDotDot() const { return _element_data->gradSlnDotDot(); }
232 
234  const FieldVariableSecond & secondSln() const { return _element_data->secondSln(Moose::Current); }
235  const FieldVariableSecond & secondSlnOld() const { return _element_data->secondSln(Moose::Old); }
237  {
238  return _element_data->secondSln(Moose::Older);
239  }
241  {
242  return _element_data->secondSln(Moose::PreviousNL);
243  }
244 
246  const FieldVariableCurl & curlSln() const { return _element_data->curlSln(Moose::Current); }
247  const FieldVariableCurl & curlSlnOld() const { return _element_data->curlSln(Moose::Old); }
248  const FieldVariableCurl & curlSlnOlder() const { return _element_data->curlSln(Moose::Older); }
249 
251  template <ComputeStage compute_stage>
253  {
254  return _element_data->template adSln<compute_stage>();
255  }
256  template <ComputeStage compute_stage>
258  {
259  return _element_data->template adGradSln<compute_stage>();
260  }
261  template <ComputeStage compute_stage>
263  {
264  return _element_data->template adSecondSln<compute_stage>();
265  }
266  template <ComputeStage compute_stage>
268  {
269  return _element_data->template adUDot<compute_stage>();
270  }
271 
273  template <ComputeStage compute_stage>
275  {
276  return _neighbor_data->template adSln<compute_stage>();
277  }
278  template <ComputeStage compute_stage>
280  {
281  return _neighbor_data->template adGradSln<compute_stage>();
282  }
283  template <ComputeStage compute_stage>
285  {
286  return _neighbor_data->template adSecondSln<compute_stage>();
287  }
288  template <ComputeStage compute_stage>
290  {
291  return _neighbor_data->template adUDot<compute_stage>();
292  }
293 
295  const FieldVariableValue & uDot() const { return _element_data->uDot(); }
296  const FieldVariableValue & uDotDot() const { return _element_data->uDotDot(); }
297  const FieldVariableValue & uDotOld() const { return _element_data->uDotOld(); }
298  const FieldVariableValue & uDotDotOld() const { return _element_data->uDotDotOld(); }
299  const VariableValue & duDotDu() const { return _element_data->duDotDu(); }
300  const VariableValue & duDotDotDu() const { return _element_data->duDotDotDu(); }
301 
304  const FieldVariableValue & slnOldNeighbor() const { return _neighbor_data->sln(Moose::Old); }
307  {
308  return _neighbor_data->sln(Moose::PreviousNL);
309  }
310 
313  {
314  return _neighbor_data->gradSln(Moose::Current);
315  }
317  {
318  return _neighbor_data->gradSln(Moose::Old);
319  }
321  {
322  return _neighbor_data->gradSln(Moose::Older);
323  }
325  {
326  return _neighbor_data->gradSln(Moose::PreviousNL);
327  }
328 
330  const FieldVariableGradient & gradSlnNeighborDot() const { return _neighbor_data->gradSlnDot(); }
332  {
333  return _neighbor_data->gradSlnDotDot();
334  }
335 
338  {
339  return _neighbor_data->secondSln(Moose::Current);
340  }
342  {
343  return _neighbor_data->secondSln(Moose::Old);
344  }
346  {
347  return _neighbor_data->secondSln(Moose::Older);
348  }
350  {
351  return _neighbor_data->secondSln(Moose::PreviousNL);
352  }
353 
356  {
357  return _neighbor_data->curlSln(Moose::Current);
358  }
360  {
361  return _neighbor_data->curlSln(Moose::Old);
362  }
364  {
365  return _neighbor_data->curlSln(Moose::Older);
366  }
367 
369  const FieldVariableValue & uDotNeighbor() const { return _neighbor_data->uDot(); }
370  const FieldVariableValue & uDotDotNeighbor() const { return _neighbor_data->uDotDot(); }
371  const FieldVariableValue & uDotOldNeighbor() const { return _neighbor_data->uDotOld(); }
372  const FieldVariableValue & uDotDotOldNeighbor() const { return _neighbor_data->uDotDotOld(); }
373  const VariableValue & duDotDuNeighbor() const { return _neighbor_data->duDotDu(); }
374  const VariableValue & duDotDotDuNeighbor() const { return _neighbor_data->duDotDotDu(); }
375 
377  template <ComputeStage compute_stage>
379  {
380  return _lower_data->template adSln<compute_stage>();
381  }
382  const FieldVariableValue & slnLower() const { return _lower_data->sln(Moose::Current); }
383 
385  virtual void computeElemValues() override;
386  virtual void computeElemValuesFace() override;
387  virtual void computeNeighborValuesFace() override;
388  virtual void computeNeighborValues() override;
389  virtual void computeLowerDValues() override;
390 
391  void setNodalValue(OutputType value, unsigned int idx = 0);
392  void setDofValues(const DenseVector<Number> & value) override;
393  Number getNodalValue(const Node & node) override;
394  Number getNodalValueOld(const Node & node) override;
395  Number getNodalValueOlder(const Node & node) override;
396  Number getElementalValue(const Elem * elem, unsigned int idx = 0) const override;
397  Number getElementalValueOld(const Elem * elem, unsigned int idx = 0) const override;
398  Number getElementalValueOlder(const Elem * elem, unsigned int idx = 0) const override;
399 
400  void getDofIndices(const Elem * elem, std::vector<dof_id_type> & dof_indices) const override;
401  const std::vector<dof_id_type> & dofIndices() const final { return _element_data->dofIndices(); }
402  unsigned int numberOfDofs() const final { return _element_data->numberOfDofs(); }
403  const std::vector<dof_id_type> & dofIndicesNeighbor() const final
404  {
405  return _neighbor_data->dofIndices();
406  }
407  const std::vector<dof_id_type> & dofIndicesLower() const final
408  {
409  return _lower_data->dofIndices();
410  }
411 
412  unsigned int numberOfDofsNeighbor() override { return _neighbor_data->dofIndices().size(); }
413 
414  void insert(NumericVector<Number> & residual) override;
415  void add(NumericVector<Number> & residual) override;
416 
417  const MooseArray<Number> & dofValue() override;
418  const MooseArray<Number> & dofValues() override;
419  const MooseArray<Number> & dofValuesOld() override;
420  const MooseArray<Number> & dofValuesOlder() override;
421  const MooseArray<Number> & dofValuesPreviousNL() override;
422  const MooseArray<Number> & dofValuesNeighbor() override;
423  const MooseArray<Number> & dofValuesOldNeighbor() override;
424  const MooseArray<Number> & dofValuesOlderNeighbor() override;
426  const MooseArray<Number> & dofValuesDot() override;
427  const MooseArray<Number> & dofValuesDotNeighbor() override;
428  const MooseArray<Number> & dofValuesDotOld() override;
429  const MooseArray<Number> & dofValuesDotOldNeighbor() override;
430  const MooseArray<Number> & dofValuesDotDot() override;
431  const MooseArray<Number> & dofValuesDotDotNeighbor() override;
432  const MooseArray<Number> & dofValuesDotDotOld() override;
434  const MooseArray<Number> & dofValuesDuDotDu() override;
435  const MooseArray<Number> & dofValuesDuDotDuNeighbor() override;
436  const MooseArray<Number> & dofValuesDuDotDotDu() override;
438 
442  template <ComputeStage compute_stage>
444 
448  void computeIncrementAtQps(const NumericVector<Number> & increment_vec);
449 
453  void computeIncrementAtNode(const NumericVector<Number> & increment_vec);
454 
461  OutputType getValue(const Elem * elem, const std::vector<std::vector<OutputType>> & phi) const;
462 
470  getGradient(const Elem * elem,
471  const std::vector<std::vector<typename OutputTools<OutputType>::OutputGradient>> &
472  grad_phi) const;
473 
477  virtual size_t phiSize() const final { return _element_data->phiSize(); }
481  virtual size_t phiFaceSize() const final { return _element_data->phiFaceSize(); }
485  virtual size_t phiNeighborSize() const final { return _neighbor_data->phiSize(); }
489  virtual size_t phiFaceNeighborSize() const final { return _neighbor_data->phiFaceSize(); }
490 
491  size_t phiLowerSize() const final { return _lower_data->phiSize(); }
492 
496  const OutputType & nodalValue();
497  const OutputType & nodalValueOld();
498  const OutputType & nodalValueOlder();
499  const OutputType & nodalValuePreviousNL();
500  const OutputType & nodalValueDot();
501  const OutputType & nodalValueDotDot();
502  const OutputType & nodalValueDotOld();
503  const OutputType & nodalValueDotDotOld();
504  const OutputType & nodalValueDuDotDu();
505  const OutputType & nodalValueDuDotDotDu();
506  const OutputType & nodalValueNeighbor();
507  const OutputType & nodalValueOldNeighbor();
508  const OutputType & nodalValueOlderNeighbor();
509  const OutputType & nodalValuePreviousNLNeighbor();
510  const OutputType & nodalValueDotNeighbor();
511  const OutputType & nodalValueDotDotNeighbor();
512  const OutputType & nodalValueDotOldNeighbor();
513  const OutputType & nodalValueDotDotOldNeighbor();
514  const OutputType & nodalValueDuDotDuNeighbor();
515  const OutputType & nodalValueDuDotDotDuNeighbor();
516 
521  {
522  return _element_data->nodalValueArray(Moose::Current);
523  }
525  {
526  return _element_data->nodalValueArray(Moose::Old);
527  }
529  {
530  return _element_data->nodalValueArray(Moose::Older);
531  }
532 
535 
536  template <ComputeStage compute_stage>
538 
539  virtual void computeNodalValues() override;
540  virtual void computeNodalNeighborValues() override;
541 
542 protected:
544 
546  std::unique_ptr<MooseVariableData<OutputType>> _element_data;
547 
549  std::unique_ptr<MooseVariableData<OutputType>> _neighbor_data;
550 
552  std::unique_ptr<MooseVariableData<OutputType>> _lower_data;
553 };
554 
555 template <typename OutputType>
556 template <ComputeStage compute_stage>
559 {
560  return _element_data->template adDofValues<compute_stage>();
561 }
562 
563 template <typename OutputType>
564 template <ComputeStage compute_stage>
567 {
568  return _element_data->template adNodalValue<compute_stage>();
569 }
const MooseArray< Number > & dofValuesDotDotOldNeighbor() override
Returns old second time derivative of neighboring degrees of freedom.
const FieldVariableValue & uDotDotNeighbor() const
Number getNodalValueOlder(const Node &node) override
Get the t-2 value of this variable at given node.
const FieldVariableCurl & curlSlnNeighbor() const
neighbor solution curls
const MooseArray< Number > & dofValuesDot() override
Returns time derivative of degrees of freedom.
const Node *const & node() const
const FieldVariableSecond & secondSlnOlderNeighbor() const
const FieldVariablePhiValue & phi() const
virtual size_t phiFaceSize() const final
Return phiFace size.
bool isNodalNeighborDefined() const
const MooseArray< Number > & dofValuesOlder() override
Returns older dof solution on element.
const FieldVariablePhiSecond & secondPhiFace() const
const FieldVariablePhiGradient & gradPhiFace() const
const OutputType & nodalValuePreviousNL()
MooseArray< std::vector< OutputDivergence > > FieldVariablePhiDivergence
const FieldVariablePhiCurl & curlPhi() const
TensorTools::IncrementRank< OutputGradient >::type OutputSecond
const FieldVariableGradient & gradSlnPreviousNL() const
MooseArray< std::vector< OutputShape > > FieldVariableTestCurl
Keeps track of stuff related to assembling.
Definition: Assembly.h:62
Class for stuff related to variables.
Definition: Adaptivity.h:29
const MooseArray< Number > & dofValuesDotOld() override
Returns old time derivative of degrees of freedom.
unsigned int TagID
Definition: MooseTypes.h:162
const FieldVariablePhiSecond & secondPhi() const
const OutputType & nodalValueOlderNeighbor()
const FieldVariableSecond & secondSlnPreviousNLNeighbor() const
const FieldVariableCurl & curlSlnOldNeighbor() const
const FieldVariablePhiSecond & secondPhiFaceNeighbor() const
const MooseArray< Number > & dofValuesDotDot() override
Returns second time derivative of degrees of freedom.
MooseArray< OutputDivergence > FieldVariableDivergence
const std::set< SubdomainID > & activeSubdomains() const override
The subdomains the variable is active on.
MooseArray< std::vector< OutputShape > > FieldVariablePhiCurl
virtual void computeLowerDValues() override
compute values at quadrature points on the lower dimensional element
const FieldVariableGradient & gradSlnNeighborDotDot() const
const FieldVariablePhiCurl & curlPhiFaceNeighbor() const
void computeIncrementAtQps(const NumericVector< Number > &increment_vec)
Compute and store incremental change in solution at QPs based on increment_vec.
bool usesPhiNeighbor() const
Whether or not this variable is actually using the shape function value.
const VariableValue & duDotDuNeighbor() const
const OutputType & nodalValueOld()
MooseArray< std::vector< OutputSecond > > FieldVariableTestSecond
const FieldVariablePhiValue & phiFace() const
MooseArray< std::vector< OutputGradient > > FieldVariablePhiGradient
void reinitAux() override
const FieldVariableSecond & secondSlnOld() const
bool activeOnSubdomain(SubdomainID subdomain) const override
Is the variable active on the subdomain?
const FieldVariableValue & slnOlderNeighbor() const
const VariableTestGradientType< OutputType, compute_stage >::type & adGradPhi()
const FieldVariablePhiSecond & secondPhiNeighbor() const
void reinitNodes(const std::vector< dof_id_type > &nodes) override
const std::vector< dof_id_type > & dofIndicesLower() const final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
const OutputType & nodalValueDotOldNeighbor()
bool isVector() const override
std::unique_ptr< MooseVariableData< OutputType > > _neighbor_data
Holder for all the data associated with the neighbor element.
void prepareNeighbor() override
Prepare the neighbor element degrees of freedom.
const OutputType & nodalValueDuDotDu()
const FieldVariableGradient & gradSlnPreviousNLNeighbor() const
const FieldVariableValue & uDotOld() const
const VariableValueType< OutputType, compute_stage >::type & adUDotNeighbor() const
const FieldVariableValue & sln() const
element solutions
const Moose::ValueType< OutputType, compute_stage >::type & adNodalValue()
virtual void prepareIC() override
Prepare the initial condition.
Base class for a system (of equations)
Definition: SystemBase.h:92
MooseArray< std::vector< OutputShape > > FieldVariableTestValue
const OutputType & nodalValuePreviousNLNeighbor()
MooseArray< OutputShape > FieldVariableValue
const FieldVariableCurl & curlSlnOld() const
MooseVariableFE(unsigned int var_num, const FEType &fe_type, SystemBase &sys, const Assembly &assembly, Moose::VarKindType var_kind, THREAD_ID tid)
const OutputType & nodalValueDotNeighbor()
const MooseArray< Number > & dofValue() override
Deprecated method.
MooseArray< OutputSecond > FieldVariableSecond
const FieldVariableGradient & gradSlnOlderNeighbor() const
virtual void computeNodalValues() override
Compute nodal values of this variable.
const FieldVariableSecond & secondSlnOldNeighbor() const
const OutputType & nodalValueDuDotDotDu()
bool usesSecondPhi() const
Whether or not this variable is computing any second derivatives.
const FieldVariablePhiCurl & curlPhiNeighbor() const
bool usesGradPhiNeighbor() const
Whether or not this variable is actually using the shape function gradient.
const FieldVariableGradient & gradSlnOldNeighbor() const
const FieldVariableValue & uDot() const
element dots
OutputType OutputShape
void prepareLowerD() override
Prepare a lower dimensional element&#39;s degrees of freedom.
const MooseArray< Number > & dofValuesDotDotOld() override
Returns old second time derivative of degrees of freedom.
const FieldVariableCurl & curlSln() const
element curls
virtual size_t phiFaceNeighborSize() const final
Return phiFaceNeighbor size.
unsigned int numberOfDofsNeighbor() override
const VariableValue & duDotDotDu() const
const FieldVariableValue & vectorTagValue(TagID tag)
const OutputType & nodalValueDot()
const FieldVariablePhiCurl & curlPhiFace() const
const FieldVariableValue & uDotDotOldNeighbor() const
const FieldVariablePhiGradient & gradPhiNeighbor() const
const MooseArray< Number > & dofValuesPreviousNL() override
Returns previous nl solution on element.
const OutputType & nodalValueDuDotDotDuNeighbor()
const FieldVariablePhiValue & phiFaceNeighbor() const
const VariableValueType< OutputType, compute_stage >::type & adSln() const
AD.
const OutputType & nodalValueNeighbor()
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
const OutputType & nodalValueDuDotDuNeighbor()
const MooseArray< Number > & dofValuesDotDotNeighbor() override
Returns second time derivative of neighboring degrees of freedom.
const FieldVariableSecond & secondSlnPreviousNL() const
void reinitAuxNeighbor() override
const MooseArray< Number > & dofValuesDotNeighbor() override
Returns time derivative of neighboring degrees of freedom.
void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const override
virtual void computeElemValuesFace() override
Compute values at facial quadrature points.
const OutputType & nodalValueDotDotNeighbor()
const FieldVariableValue & slnPreviousNL() const
unsigned int numberOfDofs() const final
Get the number of local DoFs.
const FieldVariableValue & uDotDot() const
void insert(NumericVector< Number > &residual) override
const FieldVariableSecond & secondSlnNeighbor() const
neighbor solution seconds
const FieldVariablePhiValue & phiNeighbor() const
const FieldVariableValue & slnOld() const
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:481
const MooseArray< Number > & dofValuesOld() override
Returns old dof solution on element.
const MooseArray< OutputType > & nodalValueOldArray()
const FieldVariableSecond & secondSlnOlder() const
void setDofValues(const DenseVector< Number > &value) override
Set values for this variable to keep everything up to date.
const FieldVariableValue & uDotDotOld() const
const unsigned int & currentSide() const
Current side this variable is being evaluated on.
const MooseArray< Number > & dofValuesDotOldNeighbor() override
Returns old time derivative of neighboring degrees of freedom.
virtual void computeNodalNeighborValues() override
Compute nodal values of this variable in the neighbor.
const VariableGradientType< OutputType, compute_stage >::type & adGradSlnNeighbor() const
const OutputType & nodalValue()
Methods for retrieving values of variables at the nodes.
const FieldVariableGradient & gradSlnOld() const
const Node *const & nodeNeighbor() const
TensorTools::IncrementRank< OutputShape >::type OutputGradient
size_t phiLowerSize() const final
Return the number of shape functions on the lower dimensional element for this variable.
const FieldVariableValue & uDotNeighbor() const
neighbor dots
subdomain_id_type SubdomainID
const MooseArray< OutputType > & nodalValueOlderArray()
const FieldVariableValue & increment() const
const VariableSecondType< OutputType, compute_stage >::type & adSecondSln() const
const OutputType & nodalValueOldNeighbor()
MooseArray< std::vector< OutputShape > > FieldVariablePhiValue
void reinitNodesNeighbor(const std::vector< dof_id_type > &nodes) override
std::unique_ptr< MooseVariableData< OutputType > > _lower_data
Holder for all the data associated with the lower dimeensional element.
const FieldVariableValue & slnNeighbor() const
neighbor solutions
const FieldVariableSecond & secondSln() const
element seconds
const VariableGradientType< OutputType, compute_stage >::type & adGradSln() const
const OutputType & nodalValueDotDotOldNeighbor()
void clearDofIndices() override
Clear out the dof indices.
const FieldVariableValue & slnLower() const
const MooseArray< Number > & dofValuesPreviousNLNeighbor() override
Returns previous nl solution on neighbor element.
void computeIncrementAtNode(const NumericVector< Number > &increment_vec)
Compute and store incremental change at the current node based on increment_vec.
bool isNodal() const override
Is this variable nodal.
const MooseArray< Number > & dofValuesDuDotDu() override
Returns derivative of time derivative of degrees of freedom.
virtual void computeNeighborValuesFace() override
Compute values at facial quadrature points for the neighbor.
const Elem *const & neighbor() const
Current neighboring element.
bool isNodalDefined() const
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:197
Number getElementalValueOld(const Elem *elem, unsigned int idx=0) const override
Get the old value of this variable on an element.
const VariableValueType< OutputType, compute_stage >::type & adSlnNeighbor() const
neighbor AD
MatType type
Number getElementalValueOlder(const Elem *elem, unsigned int idx=0) const override
Get the older value of this variable on an element.
const dof_id_type & nodalDofIndexNeighbor() const override
Base class for time integrators.
const dof_id_type & nodalDofIndex() const override
const FieldVariableGradient & gradSlnNeighborDot() const
neighbor grad dots
const FieldVariableGradient & gradSlnDotDot() const
OutputType getValue(const Elem *elem, const std::vector< std::vector< OutputType >> &phi) const
Compute the variable value at a point on an element.
const FieldVariableCurl & curlSlnOlder() const
const MooseArray< Number > & dofValuesNeighbor() override
Returns dof solution on neighbor element.
virtual void computeNeighborValues() override
Compute values at quadrature points for the neighbor.
const MooseArray< Number > & dofValuesDuDotDuNeighbor() override
Returns derivative of time derivative of neighboring degrees of freedom.
const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
const FieldVariableValue & matrixTagValue(TagID tag)
Number getElementalValue(const Elem *elem, unsigned int idx=0) const override
Get the current value of this variable on an element.
const OutputType & nodalValueDotOld()
void setNodalValue(OutputType value, unsigned int idx=0)
const MooseArray< Number > & dofValuesDuDotDotDu() override
Returns derivative of second time derivative of degrees of freedom.
const VariableValueType< OutputType, compute_stage >::type & adSlnLower() const
lower-d element solution
const FieldVariablePhiGradient & gradPhiLower() const
bool computingCurl() const
Whether or not this variable is computing the curl.
bool usesGradPhi() const
Whether or not this variable is actually using the shape function gradient.
TensorTools::DecrementRank< OutputShape >::type OutputDivergence
const VariableTestGradientType< OutputType, compute_stage >::type & adGradPhiFace()
const FieldVariableValue & uDotOldNeighbor() const
const FieldVariablePhiGradient & gradPhi() const
const FieldVariableValue & slnOldNeighbor() const
const OutputType & nodalValueDotDot()
const VariableValue & duDotDu() const
const MooseArray< OutputType > & nodalValueArray()
Methods for retrieving values of variables at the nodes in a MooseArray for AuxKernelBase.
MooseArray< OutputShape > FieldVariableCurl
const MooseArray< Number > & dofValuesDuDotDotDuNeighbor() override
Returns derivative of second time derivative of neighboring degrees of freedom.
MooseArray< std::vector< OutputGradient > > FieldVariableTestGradient
bool usesPhi() const
Whether or not this variable is actually using the shape function value.
const FieldVariableGradient & gradSlnOlder() const
const FieldVariableValue & slnOlder() const
void prepare() override
Prepare the elemental degrees of freedom.
const MooseArray< Number > & dofValuesOlderNeighbor() override
Returns older dof solution on neighbor element.
MooseArray< OutputGradient > FieldVariableGradient
const FieldVariableGradient & gradSlnDot() const
element gradient dots
const VariableValue & duDotDotDuNeighbor() const
void reinitNode() override
MooseArray< std::vector< OutputDivergence > > FieldVariableTestDivergence
MooseArray< std::vector< OutputSecond > > FieldVariablePhiSecond
bool usesSecondPhiNeighbor() const
Whether or not this variable is actually using the shape function second derivative on a neighbor...
const FieldVariableCurl & curlSlnOlderNeighbor() const
const MooseArray< Real > & nodalVectorTagValue(TagID tag)
const OutputType & nodalValueDotDotOld()
std::unique_ptr< MooseVariableData< OutputType > > _element_data
Holder for all the data associated with the "main" element.
Number getNodalValue(const Node &node) override
Get the value of this variable at given node.
virtual size_t phiSize() const final
Return phi size.
bool computingSecond() const
Whether or not this variable is computing any second derivatives.
TensorTools::IncrementRank< OutputShape >::type OutputGradient
Definition: MooseTypes.h:173
const MooseArray< Number > & dofValuesOldNeighbor() override
Returns old dof solution on neighbor element.
const FieldVariableGradient & gradSlnNeighbor() const
neighbor solution gradients
virtual size_t phiNeighborSize() const final
Return phiNeighbor size.
const Elem *const & currentElem() const override
Current element this variable is evaluated at.
OutputTools< OutputType >::OutputGradient getGradient(const Elem *elem, const std::vector< std::vector< typename OutputTools< OutputType >::OutputGradient >> &grad_phi) const
Compute the variable gradient value at a point on an element.
SystemBase & sys()
Get the system this variable is part of.
const OutputType & nodalValueOlder()
const VariableSecondType< OutputType, compute_stage >::type & adSecondSlnNeighbor() const
const FieldVariableGradient & gradSln() const
element gradients
OutputType OutputValue
const FieldVariableValue & slnPreviousNLNeighbor() const
virtual void computeElemValues() override
Actually compute variable values.
const Assembly & _assembly
Number getNodalValueOld(const Node &node) override
Get the old value of this variable at given node.
void add(NumericVector< Number > &residual) override
unsigned int THREAD_ID
Definition: MooseTypes.h:161
const MooseArray< Real > & nodalMatrixTagValue(TagID tag)
const FieldVariablePhiValue & phiLower() const
void prepareAux() override
const FieldVariablePhiGradient & gradPhiFaceNeighbor() const
const MooseArray< Number > & dofValues() override
Returns dof solution on element.
const MooseArray< typename Moose::RealType< compute_stage >::type > & adDofValues()
Return the AD dof values.
const VariableValueType< OutputType, compute_stage >::type & adUDot() const