https://mooseframework.inl.gov
ADMatInterfaceReaction.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 "ADMatInterfaceReaction.h"
11 
13 
16 
17 {
19  params.addParam<MaterialPropertyName>("forward_rate", "kf", "Forward reaction rate coefficient.");
20  params.addParam<MaterialPropertyName>(
21  "backward_rate", "kb", "Backward reaction rate coefficient.");
22  params.addClassDescription("Implements a reaction to establish ReactionRate=k_f*u-k_b*v "
23  "at interface.");
24  return params;
25 }
26 
28  : ADInterfaceKernel(parameters),
29  _kf(getADMaterialProperty<Real>("forward_rate")),
30  _kb(getNeighborADMaterialProperty<Real>("backward_rate"))
31 {
32 }
33 
34 ADReal
36 {
37  ADReal r = 0;
38  switch (type)
39  {
40  // Move all the terms to the LHS to get residual, for primary domain
41  // Residual = kf*u - kb*v = kf*u - kb*v
42  // Weak form for primary domain is: (test, kf*u - kb*v)
43  case Moose::Element:
44  r = _test[_i][_qp] * (_kf[_qp] * _u[_qp] - _kb[_qp] * _neighbor_value[_qp]);
45  break;
46 
47  // Similarly, weak form for secondary domain is: -(test, kf*u - kb*v),
48  // flip the sign because the direction is opposite.
49  case Moose::Neighbor:
50  r = -_test_neighbor[_i][_qp] * (_kf[_qp] * _u[_qp] - _kb[_qp] * _neighbor_value[_qp]);
51  break;
52  }
53  return r;
54 }
registerMooseObject("MooseApp", ADMatInterfaceReaction)
ADMatInterfaceReaction(const InputParameters &parameters)
const ADMaterialProperty< Real > & _kf
Forward reaction rate coefficient.
static InputParameters validParams()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
DGResidualType
Definition: MooseTypes.h:743
unsigned int _qp
Current quadrature point.
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
const ADTemplateVariableValue< T > & _u
Holds the current solution at the current quadrature point on the face.
virtual ADReal computeQpResidual(Moose::DGResidualType type) override
Compute residuals at quadrature points.
ADInterfaceKernel and ADVectorInterfaceKernel is responsible for interfacing physics across subdomain...
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
unsigned int _i
Index for test and trial functions.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
const OutputTools< T >::VariableTestValue & _test_neighbor
Side neighbor test function.
const ADTemplateVariableValue< T > & _neighbor_value
Coupled neighbor variable value.
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...
Implements a reaction to establish ReactionRate=k_f*u-k_b*v at interface.
const ADMaterialProperty< Real > & _kb
Backward reaction rate coefficient.
const OutputTools< T >::VariableTestValue & _test
Side shape function.