https://mooseframework.inl.gov
ArrayDGLowerDKernel.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 "ArrayDGKernel.h"
13 
15 
21 {
22 public:
31 
33 
37  const MooseVariableFEBase & variable() const override final { return _var; }
38 
42  const MooseVariableFEBase & lowerDVariable() const { return _lowerd_var; }
43 
47  virtual void computeResidual() override;
48 
52  virtual void computeJacobian() override;
53 
57  virtual void computeOffDiagJacobian(unsigned int jvar) override;
58 
59 protected:
63  virtual void computeLowerDResidual();
64 
69  virtual void computeLowerDQpResidual(RealEigenVector & residual) = 0;
70 
74  virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type);
75 
80 
85  const MooseVariableFEBase & jvar);
86 
91  const MooseVariableFEBase & jvar);
92 
96  virtual void initLowerDQpResidual() {}
97 
102 
108  const MooseVariableFEBase &)
109  {
110  }
111 
120 
121 private:
124 };
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type)
Computes one of the five pieces of Jacobian involving lower-d.
virtual void computeLowerDResidual()
Computes the Lower part of residual for the variable on the lower-d element.
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
RealEigenVector _work_vector
Work vector for residual computation.
This class provides an interface for common operations on field variables of both FE and FV types wit...
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
virtual void computeJacobian() override
Computes the nine pieces of element/neighbor/lower-d - element/neighbor/lower-d Jacobian.
const ArrayVariableTestValue & _test_lambda
test functions
ArrayDGLowerDKernel(const InputParameters &parameters)
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
OutputTools< RealEigenVector >::VariableValue ArrayVariableValue
Definition: MooseTypes.h:348
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
virtual void computeLowerDQpResidual(RealEigenVector &residual)=0
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
const MooseVariableFEBase & lowerDVariable() const
The variable that this kernel operates on.
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: ArrayDGKernel.h:20
forward declarations
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
virtual void computeResidual() override
Computes the residual for this element, the neighbor and the lower-d element.
virtual void computeOffDiagLowerDJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d.
virtual RealEigenVector computeLowerDQpJacobian(Moose::ConstraintJacobianType jacobian_type)=0
Computes one of the five pieces of Jacobian involving lower-d at quadrature points.
ConstraintJacobianType
Definition: MooseTypes.h:796
const ArrayMooseVariable & _lowerd_var
Variable this kernel operates on.
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
const MooseVariableFEBase & variable() const override final
The variable that this kernel operates on.
const InputParameters & parameters() const
Get the parameters of the object.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
const ArrayVariablePhiValue & _phi_lambda
Shape functions.
OutputTools< RealEigenVector >::VariableTestValue ArrayVariableTestValue
Definition: MooseTypes.h:359
const ArrayVariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
virtual RealEigenMatrix computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d at quadrature points...
ArrayMooseVariable & _var
Variable this kernel operates on.
Definition: ArrayDGKernel.h:96