https://mooseframework.inl.gov
AuxKernel.C
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 #include "AuxKernel.h"
11 
12 // local includes
13 #include "FEProblem.h"
14 #include "SubProblem.h"
15 #include "AuxiliarySystem.h"
16 #include "MooseTypes.h"
17 #include "Assembly.h"
18 #include "MortarNodalAuxKernel.h"
19 
20 #include "libmesh/numeric_vector.h"
21 #include "libmesh/dof_map.h"
22 #include "libmesh/quadrature.h"
23 #include "libmesh/boundary_info.h"
24 
25 template <typename ComputeValueType>
28 {
30 
31  if (typeid(AuxKernelTempl<ComputeValueType>).name() == typeid(VectorAuxKernel).name())
32  params.registerBase("VectorAuxKernel");
33  if (typeid(AuxKernelTempl<ComputeValueType>).name() == typeid(ArrayAuxKernel).name())
34  params.registerBase("ArrayAuxKernel");
35 
36  return params;
37 }
38 
39 template <typename ComputeValueType>
41  : AuxKernelBase(parameters),
42  MooseVariableInterface<ComputeValueType>(
43  this,
44  parameters.getCheckedPointerParam<SystemBase *>("_sys")
45  ->getVariable(parameters.get<THREAD_ID>("_tid"),
46  parameters.get<AuxVariableName>("variable"))
47  .isNodal(),
48  "variable",
50  std::is_same<Real, ComputeValueType>::value
52  : (std::is_same<RealVectorValue, ComputeValueType>::value
55 
56  _var(_aux_sys.getActualFieldVariable<ComputeValueType>(
57  _tid, parameters.get<AuxVariableName>("variable"))),
58  _nodal(_var.isNodal()),
59  _u(_nodal ? _var.nodalValueArray() : _var.sln()),
60 
61  _test(_bnd ? _var.phiFace() : _var.phi()),
62  _q_point(_bnd ? _assembly.qPointsFace() : _assembly.qPoints()),
63  _qrule(_bnd ? _assembly.qRuleFace() : _assembly.qRule()),
64  _JxW(_bnd ? _assembly.JxWFace() : _assembly.JxW()),
65  _coord(_assembly.coordTransformation()),
66 
67  _current_elem(_assembly.elem()),
68  _current_side(_assembly.side()),
69  _current_elem_volume(_assembly.elemVolume()),
70  _current_side_volume(_assembly.sideElemVolume()),
71 
72  _current_node(_assembly.node()),
73  _current_boundary_id(_assembly.currentBoundaryID()),
74  _solution(_aux_sys.solution()),
75 
76  _current_lower_d_elem(_assembly.lowerDElem())
77 {
78 }
79 
80 template <typename ComputeValueType>
81 const VariableValue &
82 AuxKernelTempl<ComputeValueType>::coupledDot(const std::string & var_name, unsigned int comp) const
83 {
84  auto var = getVar(var_name, comp);
85  if (var->kind() == Moose::VAR_AUXILIARY)
86  mooseError(
87  name(),
88  ": Unable to couple time derivative of an auxiliary variable into the auxiliary system.");
89 
90  return Coupleable::coupledDot(var_name, comp);
91 }
92 
93 template <typename ComputeValueType>
94 const VariableValue &
96  unsigned int comp) const
97 {
98  auto var = getVar(var_name, comp);
99  if (var->kind() == Moose::VAR_AUXILIARY)
100  mooseError(
101  name(),
102  ": Unable to couple time derivative of an auxiliary variable into the auxiliary system.");
103 
104  return Coupleable::coupledDotDu(var_name, comp);
105 }
106 
107 template <>
108 void
109 AuxKernelTempl<Real>::setDofValueHelper(const Real & value)
110 {
111  mooseAssert(_n_shapes == 1,
112  "Should only be calling setDofValue if there is one dof for the aux var");
113  _var.setDofValue(value, 0);
114 }
115 
116 template <>
117 void
119 {
120  mooseError("Not implemented");
121 }
122 
123 template <typename ComputeValueType>
124 void
126 {
127  if (_coincident_lower_d_calc)
128  _var.insertLower(_aux_sys.solution());
129  else
130  _var.insert(_aux_sys.solution());
131 }
132 
133 template <typename ComputeValueType>
134 void
136 {
137  precalculateValue();
138 
139  if (isNodal()) /* nodal variables */
140  {
141  mooseAssert(!_coincident_lower_d_calc,
142  "Nodal evaluations are point evaluations. We don't have to concern ourselves with "
143  "coincidence of lower-d blocks and higher-d faces because they share nodes");
144  if (_var.isNodalDefined())
145  {
146  _qp = 0;
147  ComputeValueType value = computeValue();
148  // update variable data, which is referenced by other kernels, so the value is up-to-date
149  _var.setNodalValue(value);
150  }
151  }
152  else /* elemental variables */
153  {
154  _n_shapes = _coincident_lower_d_calc ? _var.dofIndicesLower().size() : _var.numberOfDofs();
155 
156  if (_coincident_lower_d_calc)
157  {
158  static const std::string lower_error = "Make sure that the lower-d variable lives on a "
159  "lower-d block that is a superset of the boundary";
160  if (!_current_lower_d_elem)
161  mooseError("No lower-dimensional element. ", lower_error);
162  if (!_n_shapes)
163  mooseError("No degrees of freedom. ", lower_error);
164  }
165 
166  if (_n_shapes == 1) /* p0 */
167  {
168  ComputeValueType value = 0;
169  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
170  value += _JxW[_qp] * _coord[_qp] * computeValue();
171  value /= (_bnd ? _current_side_volume : _current_elem_volume);
172  if (_var.isFV())
173  setDofValueHelper(value);
174  else
175  {
176  // update the variable data referenced by other kernels.
177  // Note that this will update the values at the quadrature points too
178  // (because this is an Elemental variable)
179  if (_coincident_lower_d_calc)
180  {
181  _local_sol.resize(1);
182  if constexpr (std::is_same<Real, ComputeValueType>::value)
183  _local_sol(0) = value;
184  else
185  mooseAssert(false, "We should not enter the single dof branch with a vector variable");
186  _var.setLowerDofValues(_local_sol);
187  }
188  else
189  _var.setNodalValue(value);
190  }
191  }
192  else /* high-order */
193  {
194  _local_re.resize(_n_shapes);
195  _local_re.zero();
196  _local_ke.resize(_n_shapes, _n_shapes);
197  _local_ke.zero();
198 
199  const auto & test = _coincident_lower_d_calc ? _var.phiLower() : _test;
200 
201  // assemble the local mass matrix and the load
202  for (unsigned int i = 0; i < test.size(); i++)
203  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
204  {
205  ComputeValueType t = _JxW[_qp] * _coord[_qp] * test[i][_qp];
206  _local_re(i) += t * computeValue();
207  for (unsigned int j = 0; j < test.size(); j++)
208  _local_ke(i, j) += t * test[j][_qp];
209  }
210  // mass matrix is always SPD but in case of boundary restricted, it will be rank deficient
211  _local_sol.resize(_n_shapes);
212  if (_bnd)
213  _local_ke.svd_solve(_local_re, _local_sol);
214  else
215  _local_ke.cholesky_solve(_local_re, _local_sol);
216 
217  _coincident_lower_d_calc ? _var.setLowerDofValues(_local_sol) : _var.setDofValues(_local_sol);
218  }
219  }
220 }
221 
222 template <>
223 void
225 {
226  precalculateValue();
227 
228  if (isNodal()) /* nodal variables */
229  {
230  if (_var.isNodalDefined())
231  {
232  _qp = 0;
233  RealEigenVector value = computeValue();
234  // update variable data, which is referenced by other kernels, so the value is up-to-date
235  _var.setNodalValue(value);
236  }
237  }
238  else /* elemental variables */
239  {
240  mooseAssert(
241  _var.numberOfDofs() % _var.count() == 0,
242  "The number of degrees of freedom should be cleanly divisible by the variable count");
243  _n_shapes = _var.numberOfDofs() / _var.count();
244  if (_n_shapes == 1) /* p0 */
245  {
246  RealEigenVector value = RealEigenVector::Zero(_var.count());
247  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
248  value += _JxW[_qp] * _coord[_qp] * computeValue();
249  value /= (_bnd ? _current_side_volume : _current_elem_volume);
250  // update the variable data referenced by other kernels.
251  // Note that this will update the values at the quadrature points too
252  // (because this is an Elemental variable)
253  _var.setNodalValue(value);
254  }
255  else /* high-order */
256  {
257  _local_re.resize(_n_shapes);
258  for (unsigned int i = 0; i < _local_re.size(); ++i)
259  _local_re(i) = RealEigenVector::Zero(_var.count());
260  _local_ke.resize(_n_shapes, _n_shapes);
261  _local_ke.zero();
262 
263  // assemble the local mass matrix and the load
264  for (unsigned int i = 0; i < _test.size(); i++)
265  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
266  {
267  Real t = _JxW[_qp] * _coord[_qp] * _test[i][_qp];
268  _local_re(i) += t * computeValue();
269  for (unsigned int j = 0; j < _test.size(); j++)
270  _local_ke(i, j) += t * _test[j][_qp];
271  }
272 
273  // mass matrix is always SPD
274  _local_sol.resize(_n_shapes);
275  for (unsigned int i = 0; i < _local_re.size(); ++i)
276  _local_sol(i) = RealEigenVector::Zero(_var.count());
277  DenseVector<Number> re(_n_shapes);
278  DenseVector<Number> sol(_n_shapes);
279  for (unsigned int i = 0; i < _var.count(); ++i)
280  {
281  for (unsigned int j = 0; j < _n_shapes; ++j)
282  re(j) = _local_re(j)(i);
283 
284  if (_bnd)
285  _local_ke.svd_solve(re, sol);
286  else
287  _local_ke.cholesky_solve(re, sol);
288 
289  for (unsigned int j = 0; j < _n_shapes; ++j)
290  _local_sol(j)(i) = sol(j);
291  }
292 
293  _var.setDofValues(_local_sol);
294  }
295  }
296 }
297 
298 template <typename ComputeValueType>
301 {
302  if (_sys.solutionStatesInitialized())
303  mooseError("The solution states have already been initialized when calling ",
304  type(),
305  "::uOld().\n\n",
306  "Make sure to call uOld() within the object constructor.");
307 
308  return _nodal ? _var.nodalValueOldArray() : _var.slnOld();
309 }
310 
311 template <typename ComputeValueType>
314 {
315  if (_sys.solutionStatesInitialized())
316  mooseError("The solution states have already been initialized when calling ",
317  type(),
318  "::uOlder().\n\n",
319  "Make sure to call uOlder() within the object constructor.");
320 
321  return _nodal ? _var.nodalValueOlderArray() : _var.slnOlder();
322 }
323 
324 template <typename ComputeValueType>
325 bool
327 {
328  return dynamic_cast<MortarNodalAuxKernelTempl<ComputeValueType> *>(this) != nullptr;
329 }
330 
331 // Explicitly instantiates the three versions of the AuxKernelTempl class
332 template class AuxKernelTempl<Real>;
333 template class AuxKernelTempl<RealVectorValue>;
334 template class AuxKernelTempl<RealEigenVector>;
std::string name(const ElemQuality q)
VarFieldType
Definition: MooseTypes.h:770
Base class for auxiliary kernels.
Definition: AuxKernelBase.h:42
const OutputTools< ComputeValueType >::VariableValue & uOlder() const
Retrieves the older value of the variable that this AuxKernel operates on.
Definition: AuxKernel.C:313
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:323
virtual void compute() override
Computes the value and stores it in the solution vector.
Definition: AuxKernel.C:135
virtual const VariableValue & coupledDot(const std::string &var_name, unsigned int comp=0) const
Time derivative of a coupled variable.
Definition: Coupleable.C:1170
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Base class for a system (of equations)
Definition: SystemBase.h:84
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
unsigned int _n_shapes
number of shape functions for the finite element type and current DofObject
Definition: AuxKernel.h:148
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
bool isMortar()
Definition: AuxKernel.C:326
AuxKernelTempl(const InputParameters &parameters)
Definition: AuxKernel.C:40
virtual const OutputTools< ComputeValueType >::VariableValue & value()
The value of the variable this object is operating on.
Base class for creating new nodally-based mortar auxiliary kernels.
virtual const VariableValue & coupledDot(const std::string &var_name, unsigned int comp=0) const override
Time derivative of a coupled variable.
Definition: AuxKernel.C:82
static InputParameters validParams()
Definition: AuxKernelBase.C:20
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:343
virtual const VariableValue & coupledDotDu(const std::string &var_name, unsigned int comp=0) const
Time derivative of a coupled variable with respect to the coefficients.
Definition: Coupleable.C:1468
const OutputTools< ComputeValueType >::VariableValue & uOld() const
Retrieves the old value of the variable that this AuxKernel operates on.
Definition: AuxKernel.C:300
void insert()
Insert the just computed values into the auxiliary solution vector.
Definition: AuxKernel.C:125
MooseVariableField< ComputeValueType > & _var
This is a regular kernel so we cast to a regular MooseVariable, hides base _var.
Definition: AuxKernel.h:103
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void setDofValue(const DofValue &value, unsigned int index)=0
Degree of freedom value setters.
Interface for objects that need to get values of MooseVariables.
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
static InputParameters validParams()
Definition: AuxKernel.C:27
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
Base class for creating new auxiliary kernels and auxiliary boundary conditions.
Definition: AuxKernel.h:17
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
virtual const VariableValue & coupledDotDu(const std::string &var_name, unsigned int comp=0) const override
Time derivative of a coupled variable with respect to the coefficients.
Definition: AuxKernel.C:95
void setDofValueHelper(const ComputeValueType &dof_value)
Currently only used when the auxiliary variable is a finite volume variable, this helps call through ...
const Elem & get(const ElemType type_in)
unsigned int THREAD_ID
Definition: MooseTypes.h:237