https://mooseframework.inl.gov
MatReaction.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 "MatReaction.h"
11 
12 registerMooseObject("MooseApp", MatReaction);
13 
16 {
18  params.addCoupledVar("v",
19  "Set this to make v a coupled variable, otherwise it will use the "
20  "kernel's nonlinear variable for v");
21  params.addClassDescription("Kernel to add -L*v, where L=reaction rate, v=variable");
22  params.addRequiredParam<MaterialPropertyName>("mob_name",
23  "The reaction rate used with the kernel");
24  params.deprecateParam("mob_name", "reaction_rate", "01/01/2025");
25  params.addCoupledVar("args", "Vector of nonlinear variable arguments this object depends on");
26  return params;
27 }
28 
31  _is_coupled(isCoupled("v")),
32  _v_name(_is_coupled ? coupledName("v") : _var.name()),
33  _v(_is_coupled ? coupledValue("v") : _u),
34  _v_var(_is_coupled ? coupled("v") : _var.number()),
35  _L(getMaterialProperty<Real>("reaction_rate")),
36  _eta_name(_var.name()),
37  _dLdop(getMaterialPropertyDerivative<Real>("reaction_rate", _eta_name)),
38  _dLdv(getMaterialPropertyDerivative<Real>("reaction_rate", _v_name)),
39  _dLdarg(_n_args)
40 {
41  // Get reaction rate derivatives
42  for (unsigned int i = 0; i < _n_args; ++i)
43  _dLdarg[i] = &getMaterialPropertyDerivative<Real>("reaction_rate", i);
44 }
45 
46 void
48 {
49  validateNonlinearCoupling<Real>("reaction_rate");
50 }
51 
52 Real
54 {
55  return -_L[_qp] * _test[_i][_qp] * _v[_qp];
56 }
57 
58 Real
60 {
61  if (_is_coupled)
62  return -_dLdop[_qp] * _v[_qp] * _phi[_j][_qp] * _test[_i][_qp];
63 
64  return -(_L[_qp] + _dLdop[_qp] * _v[_qp]) * _phi[_j][_qp] * _test[_i][_qp];
65 }
66 
67 Real
69 {
70  // first handle the case where jvar is a coupled variable v being added to residual
71  // the first term in the sum just multiplies by L which is always needed
72  // the second term accounts for cases where L depends on v
73  if ((jvar == _v_var) && _is_coupled)
74  return -(_L[_qp] + _dLdv[_qp] * _v[_qp]) * _test[_i][_qp] * _phi[_j][_qp];
75 
76  // for all other vars get the coupled variable jvar is referring to
77  const unsigned int cvar = mapJvarToCvar(jvar);
78 
79  return -(*_dLdarg[cvar])[_qp] * _v[_qp] * _test[_i][_qp] * _phi[_j][_qp];
80 }
std::string name(const ElemQuality q)
static InputParameters validParams()
Definition: Kernel.C:24
static InputParameters validParams()
Definition: MatReaction.C:15
virtual void initialSetup()
Gets called at the beginning of the simulation before this object is asked to do its job...
Definition: MatReaction.C:47
This kernel adds to the residual a contribution of where is a material property and is a variable ...
Definition: MatReaction.h:22
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Interface class ("Veneer") for Kernel to provide a mapping from &#39;jvar&#39; in computeQpOffDiagJacobian in...
MatReaction(const InputParameters &parameters)
Definition: MatReaction.C:29
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
unsigned int _v_var
Definition: MatReaction.h:46
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
void deprecateParam(const std::string &old_name, const std::string &new_name, const std::string &removal_date)
const MaterialProperty< Real > & _dLdv
Reaction rate derivative w.r.t. the variable being added by this kernel.
Definition: MatReaction.h:58
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
unsigned int mapJvarToCvar(unsigned int jvar)
Return index into the _coupled_moose_vars array for a given jvar.
const MaterialProperty< Real > & _dLdop
Reaction rate derivative w.r.t. order parameter.
Definition: MatReaction.h:55
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
virtual Real computeQpResidual()
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
Definition: MatReaction.C:53
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
For coupling standard variables.
Definition: MatReaction.C:68
const bool _is_coupled
is the kernel used in a coupled form?
Definition: MatReaction.h:36
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Interface class ("Veneer") to provide generator methods for derivative material property names...
std::vector< const MaterialProperty< Real > * > _dLdarg
Reaction rate derivatives w.r.t. other coupled variables.
Definition: MatReaction.h:61
Definition: Kernel.h:15
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...
const MaterialProperty< Real > & _L
Reaction rate.
Definition: MatReaction.h:49
const VariableValue & _v
Definition: MatReaction.h:45
const unsigned int _n_args
number of coupled moose variables
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
registerMooseObject("MooseApp", MatReaction)
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: MatReaction.C:59