https://mooseframework.inl.gov
InterfaceReaction.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 "InterfaceReaction.h"
11 
13 
16 {
18  params.addRequiredParam<Real>("kf", "Forward reaction rate coefficient.");
19  params.addRequiredParam<Real>("kb", "Backward reaction rate coefficient.");
20  params.addClassDescription("Implements a reaction to establish ReactionRate=k_f*u-k_b*v "
21  "at interface.");
22  return params;
23 }
24 
26  : InterfaceKernel(parameters), _kf(getParam<Real>("kf")), _kb(getParam<Real>("kb"))
27 {
28 }
29 
30 Real
32 {
33  Real r = 0;
34  switch (type)
35  {
36  // Move all the terms to the LHS to get residual, for primary domain
37  // Residual = kf*u - kb*v = kf*u - kb*v
38  // Weak form for primary domain is: (test, kf*u - kb*v)
39  case Moose::Element:
40  r = _test[_i][_qp] * (_kf * _u[_qp] - _kb * _neighbor_value[_qp]);
41  break;
42 
43  // Similarly, weak form for secondary domain is: -(test, kf*u - kb*v),
44  // flip the sign because the direction is opposite.
45  case Moose::Neighbor:
46  r = -_test_neighbor[_i][_qp] * (_kf * _u[_qp] - _kb * _neighbor_value[_qp]);
47  break;
48  }
49  return r;
50 }
51 
52 Real
54 {
55  Real jac = 0;
56  switch (type)
57  {
59  jac = _test[_i][_qp] * _kf * _phi[_j][_qp];
60  break;
62  jac = -_test_neighbor[_i][_qp] * -_kb * _phi_neighbor[_j][_qp];
63  break;
65  jac = -_test_neighbor[_i][_qp] * _kf * _phi[_j][_qp];
66  break;
68  jac = _test[_i][_qp] * -_kb * _phi_neighbor[_j][_qp];
69  break;
70  }
71  return jac;
72 }
const TemplateVariableValue & _neighbor_value
Coupled neighbor variable value.
const TemplateVariableValue & _u
Holds the current solution at the current quadrature point on the face.
InterfaceKernel and VectorInterfaceKernel is responsible for interfacing physics across subdomains...
virtual Real computeQpJacobian(Moose::DGJacobianType type) override
Compute jacobians at quadrature points.
InterfaceReaction(const InputParameters &parameters)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const TemplateVariablePhiValue & _phi
shape function
DGResidualType
Definition: MooseTypes.h:743
static InputParameters validParams()
unsigned int _qp
Current quadrature point.
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...
Implements a reaction to establish ReactionRate=k_f*u-k_b*v at interface.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
Real _kf
Forward reaction rate coefficient.
unsigned int _i
Index for test and trial functions.
const TemplateVariableTestValue & _test
Side shape function.
virtual Real computeQpResidual(Moose::DGResidualType type) override
Compute residuals at quadrature points.
DGJacobianType
Definition: MooseTypes.h:749
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("MooseApp", InterfaceReaction)
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...
Real _kb
Backward reaction rate coefficient.
const TemplateVariableTestValue & _test_neighbor
Side neighbor test function.
const TemplateVariablePhiValue & _phi_neighbor
Side neighbor shape function.
static InputParameters validParams()