https://mooseframework.inl.gov
ADMatReaction.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 "ADMatReaction.h"
11 
13 
16 {
18  params.addClassDescription(
19  "Kernel representing the contribution of the PDE term $-L*v$, where $L$ is a reaction rate "
20  "material property, $v$ is a scalar variable (nonlinear or coupled), and whose Jacobian "
21  "contribution is calculated using automatic differentiation.");
22  params.addCoupledVar("v",
23  "Set this to make v a coupled variable, otherwise it will use the "
24  "kernel's nonlinear variable for v");
25  params.addRequiredParam<MaterialPropertyName>("reaction_rate",
26  "The reaction_rate used with the kernel.");
27  return params;
28 }
29 
31  : ADKernel(parameters),
32  _v(isCoupled("v") ? adCoupledValue("v") : _u),
33  _reaction_rate(getADMaterialProperty<Real>("reaction_rate"))
34 {
35 }
36 
37 ADReal
39 {
40  return -_reaction_rate[_qp] * _test[_i][_qp] * _v[_qp];
41 }
const ADVariableValue & _v
Kernel variable (can be nonlinear or coupled variable) (For constrained Allen-Cahn problems such as t...
Definition: ADMatReaction.h:36
virtual ADReal computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
Definition: ADMatReaction.C:38
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const ADTemplateVariableTestValue< T > & _test
the current test function
Definition: ADKernel.h:83
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
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...
const ADMaterialProperty< Real > & _reaction_rate
Reaction rate material property.
Definition: ADMatReaction.h:39
ADMatReaction(const InputParameters &parameters)
Definition: ADMatReaction.C:30
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
static InputParameters validParams()
Definition: ADMatReaction.C:15
static InputParameters validParams()
Definition: ADKernel.C:24
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
Kernel representing the contribution of the PDE term $-L*v$, where $L$ is a reaction rate coefficient...
Definition: ADMatReaction.h:19
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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...
registerMooseObject("MooseApp", ADMatReaction)
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43