www.mooseframework.org
KernelBase.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 "MooseObject.h"
13 #include "BlockRestrictable.h"
14 #include "SetupInterface.h"
16 #include "FunctionInterface.h"
17 #include "UserObjectInterface.h"
18 #include "TransientInterface.h"
19 #include "PostprocessorInterface.h"
22 #include "RandomInterface.h"
24 #include "Restartable.h"
25 #include "MeshChangedInterface.h"
26 #include "TaggingInterface.h"
27 
28 class MooseMesh;
29 class SubProblem;
30 class KernelBase;
31 class Assembly;
32 template <typename>
36 
37 template <>
39 
44 class KernelBase : public MooseObject,
45  public BlockRestrictable,
46  public SetupInterface,
48  public FunctionInterface,
49  public UserObjectInterface,
50  public TransientInterface,
54  public RandomInterface,
55  protected GeometricSearchInterface,
56  public Restartable,
57  public MeshChangedInterface,
58  public TaggingInterface
59 {
60 public:
62 
63  virtual ~KernelBase();
64 
66  virtual void computeResidual() = 0;
67 
69  virtual void computeJacobian() = 0;
70 
72  virtual void computeOffDiagJacobian(MooseVariableFEBase & jvar) = 0;
73 
74  virtual void computeADOffDiagJacobian()
75  {
76  mooseError("The computeADOffDiagJacobian method should only be called on ADKernel objects");
77  }
78 
83  virtual void computeOffDiagJacobianScalar(unsigned int jvar) = 0;
84 
89  virtual void computeNonlocalJacobian() {}
90 
95  virtual void computeNonlocalOffDiagJacobian(unsigned int /* jvar */) {}
96 
100  virtual MooseVariableFEBase & variable() = 0;
101 
106 
107 protected:
111  virtual Real computeQpJacobian() { return 0; }
116  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
117 
121  virtual void precalculateResidual() {}
122  virtual void precalculateJacobian() {}
123  virtual void precalculateOffDiagJacobian(unsigned int /* jvar */) {}
124 
125 protected:
128 
131 
134 
137 
140 
143 
144  const Elem * const & _current_elem;
145 
147  const Real & _current_elem_volume;
148 
150  unsigned int _qp;
151 
154 
156  const QBase * const & _qrule;
157 
160 
163 
165  unsigned int _i;
166 
168  unsigned int _j;
169 
172  std::vector<MooseVariableFEBase *> _save_in;
173  std::vector<AuxVariableName> _save_in_strings;
174 
177  std::vector<MooseVariableFEBase *> _diag_save_in;
178  std::vector<AuxVariableName> _diag_save_in_strings;
179 
180  std::vector<unsigned int> _displacements;
181 };
182 
Interface for objects that need parallel consistent random numbers without patterns over the course o...
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: KernelBase.h:111
virtual MooseVariableFEBase & variable()=0
Returns the variable number that this Kernel operates on.
FEProblemBase & _fe_problem
Reference to this kernel&#39;s FEProblemBase.
Definition: KernelBase.h:130
A class for creating restricted objects.
Definition: Restartable.h:29
MooseVariableFE< VectorValue< Real > > VectorMooseVariable
Definition: KernelBase.h:35
Keeps track of stuff related to assembling.
Definition: Assembly.h:62
Class for stuff related to variables.
Definition: Adaptivity.h:29
std::vector< MooseVariableFEBase * > _save_in
Definition: KernelBase.h:172
virtual void computeResidual()=0
Compute this Kernel&#39;s contribution to the residual.
virtual void precalculateOffDiagJacobian(unsigned int)
Definition: KernelBase.h:123
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:159
std::vector< AuxVariableName > _diag_save_in_strings
Definition: KernelBase.h:178
SystemBase & _sys
Reference to the EquationSystem object.
Definition: KernelBase.h:133
THREAD_ID _tid
The thread ID for this kernel.
Definition: KernelBase.h:136
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual ~KernelBase()
Definition: KernelBase.C:100
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar)=0
Computes d-residual / d-jvar... storing the result in Ke.
void mooseError(Args &&... args) const
Definition: MooseObject.h:147
Base class for a system (of equations)
Definition: SystemBase.h:92
KernelBase(const InputParameters &parameters)
Definition: KernelBase.C:59
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:176
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
virtual void precalculateJacobian()
Definition: KernelBase.h:122
MooseVariableFE< Real > MooseVariable
Definition: KernelBase.h:33
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
Definition: KernelBase.h:139
virtual void computeNonlocalOffDiagJacobian(unsigned int)
Computes d-residual / d-jvar...
Definition: KernelBase.h:95
Interface for objects that needs transient capabilities.
const Real & _current_elem_volume
Volume of the current element.
Definition: KernelBase.h:147
Interface for notifications that the mesh has changed.
This is the common base class for the two main kernel types implemented in MOOSE, EigenKernel and Ker...
Definition: KernelBase.h:44
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:42
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:65
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: KernelBase.h:177
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:156
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:74
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: KernelBase.h:171
virtual void computeNonlocalJacobian()
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries corresponding to nonlocal dofs of...
Definition: KernelBase.h:89
unsigned int _i
current index for the test function
Definition: KernelBase.h:165
Interface for objects that need to use UserObjects.
std::vector< unsigned int > _displacements
Definition: KernelBase.h:180
SubProblem & subProblem()
Returns a reference to the SubProblem for which this Kernel is active.
Definition: KernelBase.h:105
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:162
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
Definition: KernelBase.h:142
std::vector< AuxVariableName > _save_in_strings
Definition: KernelBase.h:173
unsigned int _j
current index for the shape function
Definition: KernelBase.h:168
virtual void precalculateResidual()
Following methods are used for Kernels that need to perform a per-element calculation.
Definition: KernelBase.h:121
virtual void computeJacobian()=0
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:59
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
Definition: KernelBase.h:127
An interface for accessing Materials.
Intermediate base class that ties together all the interfaces for getting MooseVariableFEBases with t...
An interface that restricts an object to subdomains via the &#39;blocks&#39; input parameter.
virtual void computeADOffDiagJacobian()
Definition: KernelBase.h:74
const Elem *const & _current_elem
Definition: KernelBase.h:144
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: KernelBase.h:116
virtual void computeOffDiagJacobianScalar(unsigned int jvar)=0
Computes jacobian block with respect to a scalar variable.
Interface for objects that need to use functions.
InputParameters validParams< KernelBase >()
Definition: KernelBase.C:22
const MooseArray< Point > & _q_point
The physical location of the element&#39;s quadrature Points, indexed by _qp.
Definition: KernelBase.h:153
Interface class for classes which interact with Postprocessors.
unsigned int THREAD_ID
Definition: MooseTypes.h:161
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:150