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.");
28 _has_coefficient(isParamValid(
"time_derivative_coefficient")),
29 _coeff(_has_coefficient && hasMaterialProperty<
Real>(
"time_derivative_coefficient")
30 ? &getMaterialProperty<
Real>(
"time_derivative_coefficient")
32 _coeff_array(_has_coefficient &&
36 _coeff_2d_array(_has_coefficient &&
43 MaterialPropertyName mat = getParam<MaterialPropertyName>(
"time_derivative_coefficient");
44 mooseError(
"Property " + mat +
" is of unsupported type for ArrayTimeDerivative");
58 "time_derivative_coefficient size is inconsistent with the number of components " 66 "time_derivative_coefficient size is inconsistent with the number of components " 69 "time_derivative_coefficient size is inconsistent with the number of components " 81 return RealEigenVector::Constant(
_var.
count(), tmp);
83 return RealEigenVector::Constant(
_var.
count(), tmp * (*_coeff)[
_qp]);
85 return tmp * (*_coeff_array)[
_qp];
87 return tmp * (*_coeff_2d_array)[
_qp].diagonal();
const MaterialProperty< Real > *const _coeff
scalar time derivative coefficient
const ArrayVariableValue & _u_dot
Time derivative of {u}.
registerMooseObject("MooseApp", ArrayTimeDerivative)
ArrayTimeDerivative(const InputParameters ¶meters)
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...
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's contribution to the diagonal Jacobian at the current quadrature point...
const ArrayVariableTestValue & _test
the current test function
unsigned int _i
current index for the test function
static InputParameters validParams()
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
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...
All array time kernels should inherit from this class.
unsigned int _j
current index for the shape function
virtual void computeQpResidual(RealEigenVector &residual) override
Compute this Kernel'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
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
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
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.
const VariableValue & _du_dot_du
Derivative of u_dot with respect to u.
unsigned int _qp
The current quadrature point index.