https://mooseframework.inl.gov
ArrayKernel.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 #include "MooseVariableScalar.h"
15 
16 class ArrayKernel : public KernelBase, public MooseVariableInterface<RealEigenVector>
17 {
18 public:
20 
22 
24  virtual void computeResidual() override;
25 
27  virtual void computeJacobian() override;
28 
30  virtual void computeOffDiagJacobian(unsigned int jvar) override;
31 
36  virtual void computeOffDiagJacobianScalar(unsigned int jvar) override;
37 
38  virtual const ArrayMooseVariable & variable() const override { return _var; }
39 
40 protected:
45  virtual void computeQpResidual(RealEigenVector & residual) = 0;
46 
51 
57 
63 
67  virtual void initQpResidual() {}
68 
72  virtual void initQpJacobian() {}
73 
78  virtual void initQpOffDiagJacobian(const MooseVariableFEBase & jvar)
79  {
80  if (jvar.number() == _var.number())
82  }
83 
86 
89 
93 
96 
99 
102 
105 
107  const unsigned int _count;
108 
109 private:
112 
115 };
virtual const ArrayMooseVariable & variable() const override
Returns the variable that this object operates on.
Definition: ArrayKernel.h:38
unsigned int number() const
Get variable number coming from libMesh.
const MappedArrayVariablePhiGradient & _array_grad_test
Definition: ArrayKernel.h:92
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void initQpOffDiagJacobian(const MooseVariableFEBase &jvar)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
Definition: ArrayKernel.h:78
This class provides an interface for common operations on field variables of both FE and FV types wit...
ArrayKernel(const InputParameters &parameters)
Definition: ArrayKernel.C:30
virtual void initQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
Definition: ArrayKernel.h:67
static InputParameters validParams()
Definition: ArrayKernel.C:23
virtual RealEigenVector computeQpJacobian()
Compute this Kernel&#39;s contribution to the diagonal Jacobian at the current quadrature point...
Definition: ArrayKernel.C:166
OutputTools< RealEigenVector >::VariableValue ArrayVariableValue
Definition: MooseTypes.h:348
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
virtual RealEigenMatrix computeQpOffDiagJacobianScalar(const MooseVariableScalar &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component...
Definition: ArrayKernel.C:241
const ArrayVariablePhiGradient & _grad_phi
gradient of the shape function
Definition: ArrayKernel.h:98
virtual void computeJacobian() override
Compute this ArrayKernel&#39;s contribution to the diagonal Jacobian entries.
Definition: ArrayKernel.C:131
std::vector< std::vector< Eigen::Map< RealDIMValue > > > MappedArrayVariablePhiGradient
Definition: MooseTypes.h:355
const ArrayVariableTestValue & _test
the current test function
Definition: ArrayKernel.h:88
const ArrayVariableTestGradient & _grad_test
gradient of the test function
Definition: ArrayKernel.h:91
virtual void computeQpResidual(RealEigenVector &residual)=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point, to be filled in residual.
const ArrayVariableGradient & _grad_u
Holds the solution gradient at current quadrature points.
Definition: ArrayKernel.h:104
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
const unsigned int _count
Number of components of the array variable.
Definition: ArrayKernel.h:107
RealEigenVector _work_vector
Work vector for residual and diag jacobian.
Definition: ArrayKernel.h:111
virtual void initQpJacobian()
Put necessary evaluations depending on qp but independent on test and shape functions here...
Definition: ArrayKernel.h:72
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component...
Definition: ArrayKernel.C:215
forward declarations
virtual void computeResidual() override
Compute this ArrayKernel&#39;s contribution to the residual.
Definition: ArrayKernel.C:95
const ArrayVariablePhiValue & _phi
the current shape functions
Definition: ArrayKernel.h:95
Interface for objects that need to get values of MooseVariables.
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ArrayKernel.C:224
Class for scalar variables (they are different).
const InputParameters & parameters() const
Get the parameters of the object.
const ArrayVariableValue & _u
Holds the solution at current quadrature points.
Definition: ArrayKernel.h:101
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
RealEigenMatrix _work_matrix
Work vector for off diag jacobian.
Definition: ArrayKernel.h:114
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes full Jacobian of jvar and the array variable this kernel operates on.
Definition: ArrayKernel.C:172
ArrayMooseVariable & _var
This is an array kernel so we cast to a ArrayMooseVariable.
Definition: ArrayKernel.h:85
OutputTools< RealEigenVector >::VariableTestGradient ArrayVariableTestGradient
Definition: MooseTypes.h:360