https://mooseframework.inl.gov
VectorFunctionReaction.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 "VectorFunctionReaction.h"
11 
13 
16 {
18  params.addClassDescription(
19  "Kernel representing the contribution of the PDE term $sign * f * u$, where $f$ "
20  "is a function coefficient and $u$ is a vector variable.");
21  MooseEnum sign("positive=1 negative=-1", "positive");
22  params.addParam<MooseEnum>("sign", sign, "Sign of boundary term in weak form.");
23  params.addParam<FunctionName>("function", "1.0", "Function coefficient multiplier for field.");
24  return params;
25 }
26 
28  : VectorKernel(parameters), _sign(getParam<MooseEnum>("sign")), _function(getFunction("function"))
29 {
30 }
31 
32 Real
34 {
35  return _sign * _function.value(_t, _q_point[_qp]) * _u[_qp] * _test[_i][_qp];
36 }
37 
38 Real
40 {
41  return _sign * _function.value(_t, _q_point[_qp]) * _phi[_j][_qp] * _test[_i][_qp];
42 }
const VectorVariableTestValue & _test
the current test function
Definition: VectorKernel.h:65
const VectorVariableValue & _u
Holds the solution at current quadrature points.
Definition: VectorKernel.h:77
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
registerMooseObject("MooseApp", VectorFunctionReaction)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
T sign(T x)
Definition: MathUtils.h:84
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
VectorFunctionReaction(const InputParameters &parameters)
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
const VectorVariablePhiValue & _phi
the current shape functions
Definition: VectorKernel.h:71
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
Definition: VectorKernel.C:21
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...
Kernel representing the contribution of the PDE term $fu$, where $f$ is a function coefficient and $u...
virtual Real value(Real t, const Point &p) const
Override this to evaluate the scalar function at point (t,x,y,z), by default this returns zero...
Definition: Function.C:44
const MooseArray< Point > & _q_point
The physical location of the element&#39;s quadrature Points, indexed by _qp.
Definition: KernelBase.h:46
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43