www.mooseframework.org
ArrayReaction.C
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 #include "ArrayReaction.h"
11 
13 
15 
18 {
20  params.addParam<MaterialPropertyName>("reaction_coefficient",
21  "The name of the reactivity, "
22  "can be scalar, vector, or matrix.");
23  params.addClassDescription("The array reaction operator with the weak "
24  "form of $(\\psi_i, u_h)$.");
25  return params;
26 }
27 
29  : ArrayKernel(parameters),
30  _r(hasMaterialProperty<Real>("reaction_coefficient")
31  ? &getMaterialProperty<Real>("reaction_coefficient")
32  : nullptr),
33  _r_array(hasMaterialProperty<RealEigenVector>("reaction_coefficient")
34  ? &getMaterialProperty<RealEigenVector>("reaction_coefficient")
35  : nullptr),
36  _r_2d_array(hasMaterialProperty<RealEigenMatrix>("reaction_coefficient")
37  ? &getMaterialProperty<RealEigenMatrix>("reaction_coefficient")
38  : nullptr)
39 {
40  if (!_r && !_r_array && !_r_2d_array)
41  {
42  MaterialPropertyName mat = getParam<MaterialPropertyName>("reaction_coefficient");
43  mooseError("Property " + mat + " is of unsupported type for ArrayReaction");
44  }
45 }
46 
49 {
50  if (_r)
51  return (*_r)[_qp] * _u[_qp] * _test[_i][_qp];
52 
53  else if (_r_array)
54  {
55  mooseAssert((*_r_array)[_qp].size() == _var.count(),
56  "reaction_coefficient size is inconsistent with the number of components of array "
57  "variable");
58  return ((*_r_array)[_qp].array() * _u[_qp].array()) * _test[_i][_qp];
59  }
60 
61  else
62  {
63  mooseAssert((*_r_2d_array)[_qp].cols() == _var.count(),
64  "reaction_coefficient size is inconsistent with the number of components of array "
65  "variable");
66  mooseAssert((*_r_2d_array)[_qp].rows() == _var.count(),
67  "reaction_coefficient size is inconsistent with the number of components of array "
68  "variable");
69  return (*_r_2d_array)[_qp] * _u[_qp] * _test[_i][_qp];
70  }
71 }
72 
75 {
76  if (_r)
77  return RealEigenVector::Constant(_var.count(), _phi[_j][_qp] * _test[_i][_qp] * (*_r)[_qp]);
78  else if (_r_array)
79  return _phi[_j][_qp] * _test[_i][_qp] * (*_r_array)[_qp];
80 
81  else
82  return _phi[_j][_qp] * _test[_i][_qp] * (*_r_2d_array)[_qp].diagonal();
83 }
84 
87 {
88  if (jvar.number() == _var.number() && _r_2d_array)
89  return _phi[_j][_qp] * _test[_i][_qp] * (*_r_2d_array)[_qp];
90  else
92 }
KernelBase::_i
unsigned int _i
current index for the test function
Definition: KernelBase.h:159
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
ArrayReaction.h
MooseObject::mooseError
void mooseError(Args &&... args) const
Definition: MooseObject.h:141
registerMooseObject
registerMooseObject("MooseApp", ArrayReaction)
ArrayReaction::validParams
static InputParameters validParams()
Definition: ArrayReaction.C:17
KernelBase::_j
unsigned int _j
current index for the shape function
Definition: KernelBase.h:162
ArrayReaction
Definition: ArrayReaction.h:19
ArrayKernel::validParams
static InputParameters validParams()
Definition: ArrayKernel.C:24
InputParameters::addParam
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object.
Definition: InputParameters.h:1198
ArrayKernel
Definition: ArrayKernel.h:21
ArrayReaction::_r_2d_array
const MaterialProperty< RealEigenMatrix > * _r_2d_array
matrix diffusion coefficient
Definition: ArrayReaction.h:36
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
KernelBase::_qp
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:144
MooseVariableBase::count
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one.
Definition: MooseVariableBase.h:98
ArrayReaction::_r_array
const MaterialProperty< RealEigenVector > * _r_array
array diffusion coefficient
Definition: ArrayReaction.h:34
ArrayReaction::ArrayReaction
ArrayReaction(const InputParameters &parameters)
Definition: ArrayReaction.C:28
InputParameters::addClassDescription
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.
Definition: InputParameters.C:70
ArrayKernel::_test
const ArrayVariableTestValue & _test
the current test function
Definition: ArrayKernel.h:102
libMesh::RealEigenMatrix
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:136
ArrayReaction::computeQpResidual
virtual RealEigenVector computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
Definition: ArrayReaction.C:48
ArrayReaction::computeQpOffDiagJacobian
virtual RealEigenMatrix computeQpOffDiagJacobian(MooseVariableFEBase &jvar) override
This is the virtual that derived classes should override for computing a full Jacobian component.
Definition: ArrayReaction.C:86
ArrayReaction::_r
const MaterialProperty< Real > * _r
scalar diffusion coefficient
Definition: ArrayReaction.h:32
ArrayReaction::computeQpJacobian
virtual RealEigenVector computeQpJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian at the current quadrature point.
Definition: ArrayReaction.C:74
ArrayKernel::computeQpOffDiagJacobian
virtual RealEigenMatrix computeQpOffDiagJacobian(MooseVariableFEBase &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component.
Definition: ArrayKernel.h:60
libMesh::RealEigenVector
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:133
ArrayKernel::_phi
const ArrayVariablePhiValue & _phi
the current shape functions
Definition: ArrayKernel.h:109
ArrayKernel::_var
ArrayMooseVariable & _var
This is an array kernel so we cast to a ArrayMooseVariable.
Definition: ArrayKernel.h:99
defineLegacyParams
defineLegacyParams(ArrayReaction)
ArrayKernel::_u
const ArrayVariableValue & _u
Holds the solution at current quadrature points.
Definition: ArrayKernel.h:115
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48