https://mooseframework.inl.gov
ArrayDGKernel.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 "DGKernelBase.h"
13 
14 #include "MooseVariableScalar.h"
15 
20 class ArrayDGKernel : public DGKernelBase, public NeighborMooseVariableInterface<RealEigenVector>
21 {
22 public:
31 
33 
37  virtual const MooseVariableFEBase & variable() const override { return _var; }
38 
42  virtual void computeOffDiagJacobian(unsigned int jvar) override;
43 
48 
53 
58  const MooseVariableFEBase & jvar) override;
59 
60 protected:
65  virtual void computeQpResidual(Moose::DGResidualType type, RealEigenVector & residual) = 0;
66 
72 
77  const MooseVariableFEBase & jvar);
78 
83 
88 
94 
126  const std::vector<Eigen::Map<RealDIMValue>> & _array_normals;
128  const unsigned int _count;
129 
130 private:
133 };
const ArrayVariablePhiGradient & _grad_phi_neighbor
Gradient of side shape function.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Override this function to consider couplings of components of the array variable. ...
const std::vector< Eigen::Map< RealDIMValue > > & _array_normals
Normals in type of Eigen map.
const ArrayVariableGradient & _grad_u_neighbor
Holds the current solution gradient at the current quadrature point.
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar) override
Computes the element-element off-diagonal Jacobian.
Class for stuff related to variables.
Definition: Adaptivity.h:31
const ArrayVariableValue & _u
Holds the current solution at the current quadrature point on the face.
Definition: ArrayDGKernel.h:98
virtual const MooseVariableFEBase & variable() const override
The variable that this kernel operates on.
Definition: ArrayDGKernel.h:37
RealEigenVector _work_vector
Work vector for residual computation.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void initQpOffDiagJacobian(Moose::DGJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
Definition: ArrayDGKernel.h:93
const ArrayVariableTestGradient & _grad_test_neighbor
Gradient of side shape function.
This class provides an interface for common operations on field variables of both FE and FV types wit...
virtual void initQpResidual(Moose::DGResidualType)
Put necessary evaluations depending on qp but independent on test functions here. ...
Definition: ArrayDGKernel.h:82
DGResidualType
Definition: MooseTypes.h:743
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override
Computes the element/neighbor-element/neighbor Jacobian.
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: ArrayDGKernel.C:27
virtual void computeElemNeighResidual(Moose::DGResidualType type) override
Computes the residual for this element or the neighbor.
const ArrayVariableTestValue & _test
test functions
virtual RealEigenMatrix computeQpOffDiagJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar)
This is the virtual that derived classes should override for computing the off-diag Jacobian...
ArrayDGKernel(const InputParameters &parameters)
Definition: ArrayDGKernel.C:33
Enhances MooseVariableInterface interface provide values from neighbor elements.
const unsigned int _count
Number of components of the array variable.
const ArrayVariableTestGradient & _grad_test
Gradient of side shape function.
virtual void initQpJacobian(Moose::DGJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
Definition: ArrayDGKernel.h:87
std::vector< std::vector< Eigen::Map< RealDIMValue > > > MappedArrayVariablePhiGradient
Definition: MooseTypes.h:355
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
const MappedArrayVariablePhiGradient & _array_grad_phi_neighbor
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: ArrayDGKernel.h:20
DGJacobianType
Definition: MooseTypes.h:749
const MappedArrayVariablePhiGradient & _array_grad_test_neighbor
forward declarations
const MappedArrayVariablePhiGradient & _array_grad_phi
virtual RealEigenVector computeQpJacobian(Moose::DGJacobianType)
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
const ArrayVariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
OutputTools< RealEigenVector >::VariablePhiGradient ArrayVariablePhiGradient
Definition: MooseTypes.h:354
virtual void computeQpResidual(Moose::DGResidualType type, RealEigenVector &residual)=0
This is the virtual that derived classes should override for computing the residual on neighboring el...
const ArrayVariablePhiValue & _phi_neighbor
Side shape function.
const InputParameters & parameters() const
Get the parameters of the object.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
OutputTools< RealEigenVector >::VariableGradient ArrayVariableGradient
Definition: MooseTypes.h:349
const MappedArrayVariablePhiGradient & _array_grad_test
const ArrayVariableTestValue & _test_neighbor
Side test function.
OutputTools< RealEigenVector >::VariableTestGradient ArrayVariableTestGradient
Definition: MooseTypes.h:360
const ArrayVariableValue & _u_neighbor
Holds the current solution at the current quadrature point.
ArrayMooseVariable & _var
Variable this kernel operates on.
Definition: ArrayDGKernel.h:96
const ArrayVariablePhiValue & _phi
Shape functions.
const ArrayVariablePhiGradient & _grad_phi
Gradient of shape function.