www.mooseframework.org
Kernel.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* DO NOT MODIFY THIS HEADER */
3 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
4 /* */
5 /* (c) 2010 Battelle Energy Alliance, LLC */
6 /* ALL RIGHTS RESERVED */
7 /* */
8 /* Prepared by Battelle Energy Alliance, LLC */
9 /* Under Contract No. DE-AC07-05ID14517 */
10 /* With the U. S. Department of Energy */
11 /* */
12 /* See COPYRIGHT for full restrictions */
13 /****************************************************************/
14 
15 #ifndef KERNEL_H
16 #define KERNEL_H
17 
18 #include "KernelBase.h"
19 
20 class Kernel;
21 
22 template <>
24 
25 class Kernel : public KernelBase
26 {
27 public:
29 
30  virtual void computeResidual() override;
31  virtual void computeJacobian() override;
32  virtual void computeOffDiagJacobian(unsigned int jvar) override;
33  virtual void computeOffDiagJacobianScalar(unsigned int jvar) override;
34 
35 protected:
37  virtual Real computeQpResidual() = 0;
38 
40  virtual Real computeQpJacobian();
41 
43  virtual Real computeQpOffDiagJacobian(unsigned int jvar);
44 
46  virtual void precalculateResidual();
47  virtual void precalculateJacobian() {}
48  virtual void precalculateOffDiagJacobian(unsigned int /* jvar */) {}
49 
51  const VariableValue & _u;
52 
55 
58 
61 };
62 
63 #endif /* KERNEL_H */
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
Definition: Kernel.C:47
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: Kernel.C:126
const VariableGradient & _grad_u
Holds the solution gradient at the current quadrature points.
Definition: Kernel.h:54
virtual Real computeQpResidual()=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
const VariableValue & _u_dot
Time derivative of u.
Definition: Kernel.h:57
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
Kernel(const InputParameters &parameters)
Definition: Kernel.C:37
const VariableValue & _du_dot_du
Derivative of u_dot with respect to u.
Definition: Kernel.h:60
This is the common base class for the two main kernel types implemented in MOOSE, EigenKernel and Ker...
Definition: KernelBase.h:47
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: Kernel.C:132
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: Kernel.C:69
InputParameters validParams< Kernel >()
Definition: Kernel.C:30
virtual void precalculateOffDiagJacobian(unsigned int)
Definition: Kernel.h:48
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:53
Definition: Kernel.h:25
virtual void precalculateResidual()
Following methods are used for Kernels that need to perform a per-element calculation.
Definition: Kernel.C:138
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: Kernel.C:97
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: Kernel.C:114
virtual void precalculateJacobian()
Definition: Kernel.h:47
const VariableValue & _u
Holds the solution at current quadrature points.
Definition: Kernel.h:51