https://mooseframework.inl.gov
VectorKernel.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 VectorKernel : public KernelBase, public MooseVariableInterface<RealVectorValue>
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 
35  virtual void computeOffDiagJacobianScalar(unsigned int jvar) override;
36 
37  virtual const VectorMooseVariable & variable() const override { return _var; }
38 
39 protected:
43  virtual Real computeQpResidual() = 0;
44 
48  virtual Real computeQpJacobian() { return 0; }
49 
54  virtual Real computeQpOffDiagJacobian(unsigned int /*jvar*/) { return 0; }
55 
59  virtual Real computeQpOffDiagJacobianScalar(unsigned int /*jvar*/) { return 0; }
60 
63 
66 
69 
72 
75 
78 
81 };
const VectorVariableTestValue & _test
the current test function
Definition: VectorKernel.h:65
virtual void computeResidual() override
Compute this VectorKernel&#39;s contribution to the residual.
Definition: VectorKernel.C:47
Class for stuff related to variables.
Definition: Adaptivity.h:31
const VectorVariableValue & _u
Holds the solution at current quadrature points.
Definition: VectorKernel.h:77
OutputTools< RealVectorValue >::VariableValue VectorVariableValue
Definition: MooseTypes.h:331
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: VectorKernel.h:54
OutputTools< RealVectorValue >::VariableTestValue VectorVariableTestValue
Definition: MooseTypes.h:341
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
OutputTools< RealVectorValue >::VariablePhiGradient VectorVariablePhiGradient
Definition: MooseTypes.h:337
const VectorVariablePhiGradient & _grad_phi
gradient of the shape function
Definition: VectorKernel.h:74
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: VectorKernel.h:48
VectorKernel(const InputParameters &parameters)
Definition: VectorKernel.C:28
virtual Real computeQpOffDiagJacobianScalar(unsigned int)
For coupling scalar variables.
Definition: VectorKernel.h:59
const VectorVariableTestGradient & _grad_test
gradient of the test function
Definition: VectorKernel.h:68
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
virtual const VectorMooseVariable & variable() const override
Returns the variable that this object operates on.
Definition: VectorKernel.h:37
virtual void computeJacobian() override
Compute this VectorKernel&#39;s contribution to the diagonal Jacobian entries.
Definition: VectorKernel.C:65
const VectorVariableGradient & _grad_u
Holds the solution gradient at current quadrature points.
Definition: VectorKernel.h:80
const VectorVariablePhiValue & _phi
the current shape functions
Definition: VectorKernel.h:71
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OutputTools< RealVectorValue >::VariableTestGradient VectorVariableTestGradient
Definition: MooseTypes.h:342
OutputTools< RealVectorValue >::VariablePhiValue VectorVariablePhiValue
Definition: MooseTypes.h:336
static InputParameters validParams()
Definition: VectorKernel.C:21
Interface for objects that need to get values of MooseVariables.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: VectorKernel.C:88
const InputParameters & parameters() const
Get the parameters of the object.
virtual Real computeQpResidual()=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
OutputTools< RealVectorValue >::VariableGradient VectorVariableGradient
Definition: MooseTypes.h:332
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: VectorKernel.C:113
const VectorMooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: VectorKernel.h:62