www.mooseframework.org
Kernel.h
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 #pragma once
11 
12 #include "KernelBase.h"
13 #include "MooseVariableInterface.h"
14 
15 class Kernel;
16 
17 template <>
19 
20 class Kernel : public KernelBase, public MooseVariableInterface<Real>
21 {
22 public:
24 
26 
28  virtual void computeResidual() override;
29 
31  virtual void computeJacobian() override;
32 
34  virtual void computeOffDiagJacobian(MooseVariableFEBase & jvar) override;
35 
39  virtual void computeOffDiagJacobian(unsigned jvar);
40 
45  virtual void computeOffDiagJacobianScalar(unsigned int jvar) override;
46 
47  virtual MooseVariable & variable() override { return _var; }
48 
49 protected:
53  virtual Real computeQpResidual() = 0;
54 
58  virtual Real computeQpJacobian() { return 0; }
59 
64  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
65 
70  {
71  return RealEigenVector::Zero(jvar.count());
72  }
73 
76 
79 
82 
85 
88 
90  const VariableValue & _u;
91 
94 };
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
Kernel::_phi
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:84
validParams< Kernel >
InputParameters validParams< Kernel >()
KernelBase.h
Kernel::Kernel
Kernel(const InputParameters &parameters)
Definition: Kernel.C:32
Kernel::computeQpOffDiagJacobianArray
virtual RealEigenVector computeQpOffDiagJacobianArray(ArrayMooseVariable &jvar)
For coupling array variables.
Definition: Kernel.h:69
MooseObject::parameters
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:76
Kernel::validParams
static InputParameters validParams()
Definition: Kernel.C:25
Kernel::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: Kernel.h:64
Kernel::variable
virtual MooseVariable & variable() override
Returns the variable that this Kernel operates on.
Definition: Kernel.h:47
VariablePhiValue
OutputTools< Real >::VariablePhiValue VariablePhiValue
Definition: MooseTypes.h:304
Kernel::computeQpJacobian
virtual Real computeQpJacobian()
Compute this Kernel's contribution to the Jacobian at the current quadrature point.
Definition: Kernel.h:58
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
Kernel::_grad_phi
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: Kernel.h:87
MooseVariableBase::count
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one.
Definition: MooseVariableBase.h:98
Kernel::computeOffDiagJacobianScalar
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: Kernel.C:203
Kernel::computeResidual
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
Definition: Kernel.C:93
VariableTestValue
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:308
Kernel::_test
const VariableTestValue & _test
the current test function
Definition: Kernel.h:78
Kernel::_grad_test
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:81
MooseArray
forward declarations
Definition: MooseArray.h:16
Kernel
Definition: Kernel.h:20
Kernel::_u
const VariableValue & _u
Holds the solution at current quadrature points.
Definition: Kernel.h:90
Kernel::computeJacobian
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
Definition: Kernel.C:113
MooseVariableInterface.h
Kernel::computeOffDiagJacobian
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: Kernel.C:135
Kernel::_var
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:75
VariableTestGradient
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:309
Kernel::computeQpResidual
virtual Real computeQpResidual()=0
Compute this Kernel's contribution to the residual at the current quadrature point.
libMesh::RealEigenVector
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:133
VariableValue
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:300
MooseVariableInterface
Interface for objects that need to get values of MooseVariables.
Definition: MooseVariableInterface.h:24
KernelBase
This is the common base class for the three main kernel types implemented in MOOSE,...
Definition: KernelBase.h:45
MooseVariableFE
Class for stuff related to variables.
Definition: Adaptivity.h:31
Kernel::_grad_u
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:93