www.mooseframework.org
MassEigenKernel.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "MassEigenKernel.h"
11 
13 
16 {
18  params.addClassDescription(
19  "An eigenkernel with weak form $\\lambda(\\psi_i, -u_h \\times coeff)$ where "
20  "$\\lambda$ is the eigenvalue.");
21  params.addParam<Real>("coefficient", 1.0, "Coefficient multiplying the term");
22  return params;
23 }
24 
26  : EigenKernel(parameters), _coef(this->template getParam<Real>("coefficient"))
27 {
28 }
29 
30 Real
32 {
33  return -_coef * _u[_qp] * _test[_i][_qp];
34 }
35 
36 Real
38 {
39  return -_coef * _phi[_j][_qp] * _test[_i][_qp];
40 }
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const Real _coef
static InputParameters validParams()
Definition: EigenKernel.C:23
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
registerMooseObject("MooseApp", MassEigenKernel)
unsigned int _i
current index for the test function
Definition: KernelBase.h:57
The behavior of this kernel is controlled by one problem-wise global parameter eigen_on_current - boo...
Definition: EigenKernel.h:23
static InputParameters validParams()
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:60
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
MassEigenKernel(const InputParameters &parameters)
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
const VariableValue & _u
Holds the solution at current quadrature points.
Definition: Kernel.h:87
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:42