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 : public KernelBase, public MooseVariableInterface<Real>
16 {
17 public:
19 
21 
23  virtual void computeResidual() override;
24 
26  virtual void computeJacobian() override;
27 
29  virtual void computeOffDiagJacobian(unsigned int jvar) override;
30 
32  virtual void computeResidualAndJacobian() override;
33 
38  virtual void computeOffDiagJacobianScalar(unsigned int jvar) override;
39 
40  virtual const MooseVariable & variable() const override { return _var; }
41 
42 protected:
46  virtual Real computeQpResidual() = 0;
47 
51  virtual Real computeQpJacobian() { return 0; }
52 
56  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
57 
61  virtual Real computeQpOffDiagJacobianScalar(unsigned int /*jvar*/) { return 0; }
62 
67  {
68  return RealEigenVector::Zero(jvar.count());
69  }
70 
73 
76 
79 
82 
85 
87  const VariableValue & _u;
88 
91 };
OutputTools< Real >::VariableGradient VariableGradient
Definition: MooseTypes.h:303
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
Definition: Kernel.C:91
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:90
static InputParameters validParams()
Definition: Kernel.C:23
virtual Real computeQpResidual()=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
virtual RealEigenVector computeQpOffDiagJacobianArray(const ArrayMooseVariable &jvar)
For coupling array variables.
Definition: Kernel.h:66
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:72
Class for stuff related to variables.
Definition: Adaptivity.h:31
virtual Real computeQpOffDiagJacobian(unsigned int)
For coupling standard variables.
Definition: Kernel.h:56
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: Kernel.h:84
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeResidualAndJacobian() override
Compute the residual and Jacobian together.
Definition: Kernel.C:198
Kernel(const InputParameters &parameters)
Definition: Kernel.C:30
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:22
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:312
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: Kernel.h:51
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: Kernel.C:114
OutputTools< Real >::VariablePhiGradient VariablePhiGradient
Definition: MooseTypes.h:308
forward declarations
Definition: MooseArray.h:17
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual const MooseVariable & variable() const override
Returns the variable that this object operates on.
Definition: Kernel.h:40
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:78
Interface for objects that need to get values of MooseVariables.
Definition: Kernel.h:15
const InputParameters & parameters() const
Get the parameters of the object.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:138
virtual Real computeQpOffDiagJacobianScalar(unsigned int)
For coupling scalar variables.
Definition: Kernel.h:61
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:313
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: Kernel.C:136
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: Kernel.C:184
const VariableValue & _u
Holds the solution at current quadrature points.
Definition: Kernel.h:87