https://mooseframework.inl.gov
DGLowerDKernel.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 "DGKernel.h"
17 class DGLowerDKernel : public DGKernel
18 {
19 public:
28 
30 
34  const MooseVariableFEBase & variable() const override final { return _var; }
35 
39  const MooseVariableFEBase & lowerDVariable() const { return _lowerd_var; }
40 
44  virtual void computeResidual() override;
45 
49  virtual void computeJacobian() override;
50 
54  virtual void computeOffDiagJacobian(unsigned int jvar) override;
55 
56 protected:
60  virtual void computeLowerDResidual();
61 
66  virtual Real computeLowerDQpResidual() = 0;
67 
71  virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type);
72 
77 
82  const MooseVariableFEBase & jvar);
83 
88  const MooseVariableFEBase & /*jvar*/)
89  {
90  return 0;
91  }
92 
96  virtual void initLowerDQpResidual() {}
97 
102 
108  const MooseVariableFEBase &)
109  {
110  }
111 
120 };
const MooseVariable & _lowerd_var
Variable this kernel operates on.
virtual void computeOffDiagLowerDJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d.
virtual Real computeLowerDQpResidual()=0
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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.
const MooseVariableFEBase & variable() const override final
The variable that this kernel operates on.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
MooseVariable & _var
Variable this kernel operates on.
Definition: DGKernel.h:92
const VariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
const VariablePhiValue & _phi_lambda
Shape functions.
const VariableTestValue & _test_lambda
test functions
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:18
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:324
virtual void computeResidual() override
Computes the residual for this element, the neighbor and the lower-d element.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
virtual void computeJacobian() override
Computes the nine pieces of element/neighbor/lower-d - element/neighbor/lower-d Jacobian.
virtual void computeLowerDResidual()
Computes the Lower part of residual for the variable on the lower-d element.
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType jacobian_type)=0
Computes one of the five pieces of Jacobian involving lower-d at quadrature points.
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type)
Computes one of the five pieces of Jacobian involving lower-d.
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:314
forward declarations
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
ConstraintJacobianType
Definition: MooseTypes.h:796
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
virtual Real computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d at quadrature points...
DGLowerDKernel(const InputParameters &parameters)
const InputParameters & parameters() const
Get the parameters of the object.
const MooseVariableFEBase & lowerDVariable() const
The variable that this kernel operates on.
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...