https://mooseframework.inl.gov
ArrayTimeDerivative.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 "ArrayTimeDerivative.h"
11 
13 
16 {
18  params.addClassDescription("Array time derivative operator with the weak form of $(\\psi_i, "
19  "\\frac{\\partial u_h}{\\partial t})$.");
20  params.addParam<MaterialPropertyName>("time_derivative_coefficient",
21  "The name of the time derivative coefficient. "
22  "Can be scalar, vector, or matrix material property.");
23  return params;
24 }
25 
27  : ArrayTimeKernel(parameters),
28  _has_coefficient(isParamValid("time_derivative_coefficient")),
29  _coeff(_has_coefficient && hasMaterialProperty<Real>("time_derivative_coefficient")
30  ? &getMaterialProperty<Real>("time_derivative_coefficient")
31  : nullptr),
32  _coeff_array(_has_coefficient &&
33  hasMaterialProperty<RealEigenVector>("time_derivative_coefficient")
34  ? &getMaterialProperty<RealEigenVector>("time_derivative_coefficient")
35  : nullptr),
36  _coeff_2d_array(_has_coefficient &&
37  hasMaterialProperty<RealEigenMatrix>("time_derivative_coefficient")
38  ? &getMaterialProperty<RealEigenMatrix>("time_derivative_coefficient")
39  : nullptr)
40 {
42  {
43  MaterialPropertyName mat = getParam<MaterialPropertyName>("time_derivative_coefficient");
44  mooseError("Property " + mat + " is of unsupported type for ArrayTimeDerivative");
45  }
46 }
47 
48 void
50 {
51  if (!_has_coefficient)
52  residual = _u_dot[_qp] * _test[_i][_qp];
53  else if (_coeff)
54  residual = (*_coeff)[_qp] * _u_dot[_qp] * _test[_i][_qp];
55  else if (_coeff_array)
56  {
57  mooseAssert((*_coeff_array)[_qp].size() == _var.count(),
58  "time_derivative_coefficient size is inconsistent with the number of components "
59  "in array variable");
60  // WARNING: use noalias() syntax with caution. See ArrayDiffusion.C for more details.
61  residual.noalias() = (*_coeff_array)[_qp].asDiagonal() * _u_dot[_qp] * _test[_i][_qp];
62  }
63  else
64  {
65  mooseAssert((*_coeff_2d_array)[_qp].cols() == _var.count(),
66  "time_derivative_coefficient size is inconsistent with the number of components "
67  "in array variable");
68  mooseAssert((*_coeff_2d_array)[_qp].rows() == _var.count(),
69  "time_derivative_coefficient size is inconsistent with the number of components "
70  "in array variable");
71  // WARNING: use noalias() syntax with caution. See ArrayDiffusion.C for more details.
72  residual.noalias() = (*_coeff_2d_array)[_qp] * _u_dot[_qp] * _test[_i][_qp];
73  }
74 }
75 
78 {
79  Real tmp = _test[_i][_qp] * _phi[_j][_qp] * _du_dot_du[_qp];
80  if (!_has_coefficient)
81  return RealEigenVector::Constant(_var.count(), tmp);
82  else if (_coeff)
83  return RealEigenVector::Constant(_var.count(), tmp * (*_coeff)[_qp]);
84  else if (_coeff_array)
85  return tmp * (*_coeff_array)[_qp];
86  else
87  return tmp * (*_coeff_2d_array)[_qp].diagonal();
88 }
89 
92 {
93  if (jvar.number() == _var.number() && _coeff_2d_array)
94  return _phi[_j][_qp] * _test[_i][_qp] * _du_dot_du[_qp] * (*_coeff_2d_array)[_qp];
95  else
97 }
const MaterialProperty< Real > *const _coeff
scalar time derivative coefficient
const ArrayVariableValue & _u_dot
Time derivative of {u}.
registerMooseObject("MooseApp", ArrayTimeDerivative)
ArrayTimeDerivative(const InputParameters &parameters)
unsigned int number() const
Get variable number coming from libMesh.
static InputParameters validParams()
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const MaterialProperty< RealEigenMatrix > *const _coeff_2d_array
matrix time derivative coefficient
This class provides an interface for common operations on field variables of both FE and FV types wit...
virtual RealEigenVector computeQpJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian at the current quadrature point...
const ArrayVariableTestValue & _test
the current test function
Definition: ArrayKernel.h:88
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
static InputParameters validParams()
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
const MaterialProperty< RealEigenVector > *const _coeff_array
array time derivative coefficient
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component...
Definition: ArrayKernel.C:215
All array time kernels should inherit from this class.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
virtual void computeQpResidual(RealEigenVector &residual) override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point, to be filled in residual.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const ArrayVariablePhiValue & _phi
the current shape functions
Definition: ArrayKernel.h:95
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase &jvar) override
This is the virtual that derived classes should override for computing a full Jacobian component...
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
const bool _has_coefficient
whether or not the coefficient property is provided
ArrayMooseVariable & _var
This is an array kernel so we cast to a ArrayMooseVariable.
Definition: ArrayKernel.h:85
const VariableValue & _du_dot_du
Derivative of u_dot with respect to u.
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43