https://mooseframework.inl.gov
MooseVariableDataFV.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 "MooseArray.h"
13 #include "MooseTypes.h"
14 #include "MeshChangedInterface.h"
15 #include "MooseVariableDataBase.h"
16 #include "TheWarehouse.h"
17 
18 #include "libmesh/tensor_tools.h"
19 #include "libmesh/vector_value.h"
20 #include "libmesh/tensor_value.h"
21 #include "libmesh/type_n_tensor.h"
22 #include "libmesh/fe_type.h"
23 #include "libmesh/dof_map.h"
24 #include "libmesh/enum_fe_family.h"
25 #include "SubProblem.h"
26 
27 #include <functional>
28 #include <vector>
29 
30 class FaceInfo;
31 class SystemBase;
32 class TimeIntegrator;
33 class Assembly;
34 
35 template <typename>
36 class MooseVariableFV;
37 
38 namespace libMesh
39 {
40 class QBase;
41 }
42 
43 namespace Moose
44 {
45 template <typename T>
46 void
47 initDofIndices(T & data, const Elem & elem)
48 {
49  if (data._prev_elem != &elem)
50  {
51  data._dof_map.dof_indices(&elem, data._dof_indices, data._var_num);
52  data._prev_elem = &elem;
53  }
54 }
55 }
56 
57 namespace
58 {
59 template <typename T, typename T2>
60 void
61 assignForAllQps(const T & value, T2 & array, const unsigned int nqp)
62 {
63  for (const auto qp : make_range(nqp))
64  array[qp] = value;
65 }
66 }
67 
68 template <typename OutputType>
70 {
71 public:
72  // type for gradient, second and divergence of template class OutputType
76 
77  // shortcut for types storing values on quadrature points
83 
84  // shape function type for the template class OutputType
86 
87  // type for gradient, second and divergence of shape functions of template class OutputType
91 
92  // DoF value type for the template class OutputType
95 
97  SystemBase & sys,
98  THREAD_ID tid,
99  Moose::ElementType element_type,
100  const Elem * const & elem);
101 
102  bool isNodal() const override { return false; }
103  bool hasDoFsOnNodes() const override { return false; }
105 
109  bool needsAD() const { return _need_ad; }
110 
115  void setGeometry(Moose::GeometryType gm_type);
116 
118 
122  void computeValuesFace(const FaceInfo & fi);
123 
127  void computeValues();
128 
132  void computeAD(const unsigned int num_dofs, const unsigned int nqp);
133 
135 
139  const Elem * const & currentElem() const { return _elem; }
140 
144  void prepareIC();
145 
147 
151  const FieldVariableValue & sln(Moose::SolutionState state) const;
152 
156  const FieldVariableGradient & gradSlnDot() const;
157 
161  const FieldVariableGradient & gradSlnDotDot() const;
162 
168 
173  const FieldVariableCurl & curlSln(Moose::SolutionState state) const;
174 
176  {
177  _need_ad = _need_ad_u = true;
178  return _ad_u;
179  }
180 
182  {
183  _need_ad = _need_ad_grad_u = true;
184  return _ad_grad_u;
185  }
186 
188  {
189  mooseError("Gradient of time derivative not yet implemented for FV");
190  }
191 
193  {
194  _need_ad = _need_ad_second_u = true;
195  return _ad_second_u;
196  }
197 
199 
201 
202  const FieldVariableValue & uDot() const;
203 
204  const FieldVariableValue & uDotDot() const;
205 
206  const FieldVariableValue & uDotOld() const;
207 
208  const FieldVariableValue & uDotDotOld() const;
209 
210  const VariableValue & duDotDu() const
211  {
212  _need_du_dot_du = true;
213  return _du_dot_du;
214  }
215 
216  const VariableValue & duDotDotDu() const
217  {
218  _need_du_dotdot_du = true;
219  return _du_dotdot_du;
220  }
221 
225  void setDofValues(const DenseVector<OutputData> & values);
226 
228 
231  void setDofValue(const OutputData & value, unsigned int index);
233 
234  OutputData
235  getElementalValue(const Elem * elem, Moose::SolutionState state, unsigned int idx = 0) const;
236 
238 
239  void getDofIndices(const Elem * elem, std::vector<dof_id_type> & dof_indices) const;
240  const std::vector<dof_id_type> & dofIndices() const;
241  unsigned int numberOfDofs() const;
243  {
244  _dof_indices.clear();
245  _prev_elem = nullptr;
246  }
247 
249 
250  const DoFValue & dofValuesDot() const;
251  const DoFValue & dofValuesDotOld() const;
252  const DoFValue & dofValuesDotDot() const;
253  const DoFValue & dofValuesDotDotOld() const;
256 
260  const MooseArray<ADReal> & adDofValues() const;
261 
265  const MooseArray<ADReal> & adDofValuesDot() const;
266 
268 
273  const FieldVariableValue & increment() const { return _increment; }
274 
279 
281  bool hasDirichletBC() const { return _has_dirichlet_bc; }
282 
283  void meshChanged() override;
284 
285 protected:
286  virtual const MooseVariableFV<OutputType> & var() const override { return _var; }
287 
288 private:
289  void initializeSolnVars();
290 
296  void fetchADDoFValues();
297 
301  bool safeToComputeADUDot() const;
302 
305 
307 
308  const unsigned int _var_num;
309 
311 
314 
317 
320 
323 
325  mutable bool _need_second;
326  mutable bool _need_second_old;
327  mutable bool _need_second_older;
329 
331  mutable bool _need_curl;
332  mutable bool _need_curl_old;
333  mutable bool _need_curl_older;
334 
336  mutable bool _need_ad;
337  mutable bool _need_ad_u;
338  mutable bool _need_ad_u_dot;
339  mutable bool _need_ad_u_dotdot;
340  mutable bool _need_ad_grad_u;
341  mutable bool _need_ad_grad_u_dot;
342  mutable bool _need_ad_second_u;
343 
347 
353 
358 
369 
370  // time derivatives
371 
374 
377 
380 
383 
386 
389 
392 
396  const Elem * const & _elem;
398  mutable const Elem * _prev_elem = nullptr;
399 
400  const std::vector<dof_id_type> & initDofIndices();
401 
404 
406  const bool _displaced;
407 
410 
413 
418 
467 
468  friend void Moose::initDofIndices<>(MooseVariableDataFV<OutputType> &, const Elem &);
469 };
470 
472 
473 template <typename OutputType>
474 const MooseArray<ADReal> &
476 {
477  _need_ad = true;
478  return _ad_dof_values;
479 }
480 
481 template <typename OutputType>
482 const MooseArray<ADReal> &
484 {
485  _need_ad = _need_ad_u_dot = true;
486  return _ad_dofs_dot;
487 }
488 
489 template <typename OutputType>
490 inline bool
492 {
493  // If we don't have a time integrator then we have no way to calculate _ad_u_dot because we rely
494  // on calls to TimeIntegrator::computeADTimeDerivatives. Another potential situation where
495  // _ad_u_dot computation is potentially troublesome is if we are an auxiliary variable which uses
496  // the auxiliary system copy of the time integrator. Some derived time integrator classes do setup
497  // in their solve() method, and that solve() method only happens for the nonlinear system copy of
498  // the time integrator.
499  return _time_integrator && (_var.kind() == Moose::VAR_SOLVER);
500 }
501 
502 template <typename OutputType>
505 {
506  _need_ad = _need_ad_u_dot = true;
507 
508  if (!safeToComputeADUDot())
509  // We will just copy the value of _u_dot into _ad_u_dot
510  _need_u_dot = true;
511 
512  return _ad_u_dot;
513 }
514 
515 template <typename OutputType>
518 {
519  // Generally speaking, we need u dot information when computing u dot dot
520  adUDot();
521 
522  _need_ad = _need_ad_u_dotdot = true;
523 
524  if (!safeToComputeADUDot())
525  // We will just copy the value of _u_dotdot into _ad_u_dotdot
526  _need_u_dotdot = true;
527 
528  return _ad_u_dotdot;
529 }
530 
531 template <typename OutputType>
532 const std::vector<dof_id_type> &
534 {
535  return const_cast<MooseVariableDataFV<OutputType> *>(this)->initDofIndices();
536 }
537 
538 template <typename OutputType>
539 unsigned int
541 {
542  return dofIndices().size();
543 }
void initDofIndices(T &data, const Elem &elem)
typename OutputTools< typename Moose::ADType< T >::type >::VariableSecond ADTemplateVariableSecond
Definition: MooseTypes.h:609
const VariableValue & duDotDu() const
libMesh::FEContinuity getContinuity() const override
Return the variable continuity.
const unsigned int _var_num
FieldVariableValue _u_dot
u_dot (time derivative)
const std::vector< dof_id_type > & initDofIndices()
FieldVariableGradient _grad_u_dotdot
const FieldVariableValue & uDotOld() const
VariableValue _du_dotdot_du
derivative of u_dotdot wrt u
QueryCache is a convenient way to construct and pass around (possible partially constructed) warehous...
Definition: TheWarehouse.h:208
FieldVariableValue _u_dotdot_old
u_dotdot_old (second time derivative)
const MooseArray< ADReal > & adDofValuesDot() const
Return the AD dof time derivatives.
Keeps track of stuff related to assembling.
Definition: Assembly.h:101
void computeIncrementAtQps(const libMesh::NumericVector< libMesh::Number > &increment_vec)
Compute and store incremental change in solution at QPs based on increment_vec.
bool needsAD() const
Returns whether this data structure needs automatic differentiation calculations. ...
const ADTemplateVariableSecond< OutputType > & adSecondSln() const
bool hasDoFsOnNodes() const override
Whether this data is associated with a variable that has DoFs on nodes.
const FieldVariableValue & uDot() const
const MooseVariableFV< OutputType > & _var
A const reference to the owning MooseVariableFV object.
const MooseArray< ADReal > & adDofValues() const
Return the AD dof values.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
void setDofValue(const OutputData &value, unsigned int index)
dof value setters
OutputType type
Definition: MooseTypes.h:268
OutputType type
Definition: MooseTypes.h:257
bool hasDirichletBC() const
checks if a Dirichlet BC exists on this face
const bool _displaced
Whether this variable is being calculated on a displaced system.
unsigned int numberOfDofs() const
const MooseArray< libMesh::Number > & dofValuesDuDotDu() const
const FieldVariableValue & uDotDotOld() const
const ADTemplateVariableValue< OutputType > & adSln() const
void setGeometry(Moose::GeometryType gm_type)
Set the geometry type before calculating variables values.
const Assembly & _assembly
MooseArray< ADReal > _ad_dof_values
virtual const MooseVariableFV< OutputType > & var() const override
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:603
bool _has_dirichlet_bc
if this variable has a dirichlet bc defined on a particular face
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
MooseArray< OutputSecond > FieldVariableSecond
ElementType
Definition: MooseTypes.h:763
Base class for a system (of equations)
Definition: SystemBase.h:84
void fetchADDoFValues()
Helper methods for assigning nodal values from their corresponding solution values (dof values as the...
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
A structure for storing the various lists that contain the names of the items to be exported...
FieldVariableSecond _second_u
second_u
FieldVariableCurl _curl_u
curl_u
const TimeIntegrator *const _time_integrator
Pointer to time integrator.
TheWarehouse::QueryCache _fv_flux_kernel_query_cache
Cached warehouse query for FVFluxKernels.
void getDofIndices(const Elem *elem, std::vector< dof_id_type > &dof_indices) const
MooseArray< ADReal > _ad_dofs_dot
const ADTemplateVariableGradient< OutputType > & adGradSln() const
const ADTemplateVariableValue< OutputType > & adUDot() const
FieldVariableCurl _curl_u_older
const ADTemplateVariableGradient< OutputType > & adGradSlnDot() const
const DoFValue & dofValuesDotDot() const
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
libMesh::TensorTools::IncrementRank< OutputShape >::type OutputShapeGradient
FieldVariableValue _u_dotdot
u_dotdot (second time derivative)
const FieldVariableValue & sln(Moose::SolutionState state) const
Local solution value.
Interface for notifications that the mesh has changed.
ADTemplateVariableGradient< OutputShape > _ad_grad_u_dot
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
ADTemplateVariableValue< OutputShape > _ad_u_dot
libMesh::TensorTools::IncrementRank< OutputGradient >::type OutputSecond
void setDofValues(const DenseVector< OutputData > &values)
Set local DOF values and evaluate the values on quadrature points.
const VariableValue & duDotDotDu() const
MooseArray< OutputType > FieldVariableCurl
const libMesh::FEType & _fe_type
void prepareIC()
prepare the initial condition
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:606
ADReal _ad_real_dummy
A dummy ADReal variable.
const FieldVariableValue & uDotDot() const
const DoFValue & dofValuesDot() const
FieldVariableSecond _second_u_older
const DoFValue & dofValuesDotOld() const
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu() const
const FieldVariableValue & increment() const
Increment getter.
bool isNodal() const override
const DoFValue & dofValuesDotDotOld() const
void computeValues()
compute the variable values
FieldVariableSecond _second_u_old
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:314
libMesh::TensorTools::IncrementRank< OutputShapeGradient >::type OutputShapeSecond
FieldVariableValue _increment
Increment in the variable used in dampers.
FieldVariableSecond _second_u_previous_nl
Moose::DOFType< OutputType >::type OutputData
libMesh::TensorTools::IncrementRank< OutputType >::type OutputGradient
const libMesh::QBase * _qrule
The quadrature rule.
Base class for time integrators.
Moose::ElementType _element_type
The element type this object is storing data for. This is either Element, Neighbor, or Lower.
GeometryType
Definition: MooseTypes.h:248
void meshChanged() override
Called on this object when the mesh changes.
Moose::ShapeType< OutputType >::type OutputShape
libMesh::TensorTools::DecrementRank< OutputType >::type OutputDivergence
SolutionState
Definition: MooseTypes.h:233
bool _need_second
SolutionState second_u flags.
FieldVariableGradient _grad_u_dot
grad_u dots
libMesh::FEContinuity _continuity
Continuity type of the variable.
void computeValuesFace(const FaceInfo &fi)
compute the variable values
MooseVariableDataFV(const MooseVariableFV< OutputType > &var, SystemBase &sys, THREAD_ID tid, Moose::ElementType element_type, const Elem *const &elem)
bool safeToComputeADUDot() const
Helper method that tells us whether it&#39;s safe to compute _ad_u_dot.
IntRange< T > make_range(T beg, T end)
const std::vector< dof_id_type > & dofIndices() const
OutputData getElementalValue(const Elem *elem, Moose::SolutionState state, unsigned int idx=0) const
libMesh::TensorTools::DecrementRank< OutputShape >::type OutputShapeDivergence
MooseArray< OutputGradient > FieldVariableGradient
FieldVariableValue _u_dot_old
u_dot_old (time derivative)
ADTemplateVariableSecond< OutputShape > _ad_second_u
MooseArray< OutputDivergence > FieldVariableDivergence
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
MooseArray< ADReal > _ad_dofs_dotdot
const Elem *const & _elem
The current elem.
MooseArray< OutputData > DoFValue
ADTemplateVariableValue< OutputShape > _ad_u
AD u.
This class provides variable solution values for other classes/objects to bind to when looping over f...
const FieldVariableCurl & curlSln(Moose::SolutionState state) const
Local solution curl getter.
const ADTemplateVariableValue< OutputType > & adUDotDot() const
VariableValue _du_dot_du
derivative of u_dot wrt u
const FieldVariableSecond & secondSln(Moose::SolutionState state) const
Local solution second spatial derivative getter.
FieldVariableCurl _curl_u_old
std::vector< dof_id_type > _dof_indices
The dof indices for the current element.
void computeAD(const unsigned int num_dofs, const unsigned int nqp)
compute AD things
const FieldVariableGradient & gradSlnDot() const
Local time derivative of solution gradient getter.
MooseArray< OutputType > FieldVariableValue
ADTemplateVariableGradient< OutputShape > _ad_grad_u
bool _need_curl
curl flags
const FieldVariableGradient & gradSlnDotDot() const
Local second time derivative of solution gradient getter.
const Elem * _prev_elem
used to keep track of when dof indices are out of date
unsigned int THREAD_ID
Definition: MooseTypes.h:209
TheWarehouse::QueryCache _fv_elemental_kernel_query_cache
Cached warehouse query for FVElementalKernels.
ADTemplateVariableValue< OutputShape > _ad_u_dotdot
const Elem *const & currentElem() const
The current element.
const ADReal _ad_zero
A zero AD variable.