www.mooseframework.org
DGKernel.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
19 {
20 public:
29 
31 
35  virtual const MooseVariableFEBase & variable() const override { return _var; }
36 
41 
46 
51  const MooseVariableFEBase & jvar) override;
52 
53 protected:
59 
65 
69  virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar);
70 
75 
81 
87  const MooseVariableFEBase & /*jvar*/)
88  {
89  }
90 
94  const VariableValue & _u;
117 };
const VariablePhiGradient & _grad_phi_neighbor
Gradient of side shape function.
Definition: DGKernel.h:108
OutputTools< Real >::VariableGradient VariableGradient
Definition: MooseTypes.h:303
const VariableTestGradient & _grad_test
Gradient of side shape function.
Definition: DGKernel.h:104
virtual const MooseVariableFEBase & variable() const override
The variable that this kernel operates on.
Definition: DGKernel.h:35
Class for stuff related to variables.
Definition: Adaptivity.h:31
const VariableValue & _u_neighbor
Holds the current solution at the current quadrature point.
Definition: DGKernel.h:114
virtual void computeElemNeighResidual(Moose::DGResidualType type) override
Computes the residual for this element or the neighbor.
Definition: DGKernel.C:111
const VariablePhiValue & _phi_neighbor
Side shape function.
Definition: DGKernel.h:106
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real computeQpJacobian(Moose::DGJacobianType type)=0
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
const VariableTestGradient & _grad_test_neighbor
Gradient of side shape function.
Definition: DGKernel.h:112
This class provides an interface for common operations on field variables of both FE and FV types wit...
DGResidualType
Definition: MooseTypes.h:656
const VariablePhiGradient & _grad_phi
Gradient of shape function.
Definition: DGKernel.h:100
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
MooseVariable & _var
Variable this kernel operates on.
Definition: DGKernel.h:92
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernel.C:27
const VariableTestValue & _test_neighbor
Side test function.
Definition: DGKernel.h:110
Enhances MooseVariableInterface interface provide values from neighbor elements.
const VariableGradient & _grad_u_neighbor
Holds the current solution gradient at the current quadrature point.
Definition: DGKernel.h:116
virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar)
This is the virtual that derived classes should override for computing the off-diag Jacobian...
Definition: DGKernel.C:218
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:18
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:312
virtual void precalculateQpOffDiagJacobian(Moose::DGJacobianType, const MooseVariableFEBase &)
Insertion point for evaluations that depend on qp but are independent of the test and shape functions...
Definition: DGKernel.h:86
DGKernel(const InputParameters &parameters)
Definition: DGKernel.C:33
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
virtual void precalculateQpJacobian(Moose::DGJacobianType)
Insertion point for evaluations that depend on qp but are independent of the test and shape functions...
Definition: DGKernel.h:80
virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override
Computes the element/neighbor-element/neighbor Jacobian.
Definition: DGKernel.C:148
const VariableTestValue & _test
test functions
Definition: DGKernel.h:102
DGJacobianType
Definition: MooseTypes.h:662
OutputTools< Real >::VariablePhiGradient VariablePhiGradient
Definition: MooseTypes.h:308
forward declarations
Definition: MooseArray.h:17
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
Definition: DGKernel.h:96
const InputParameters & parameters() const
Get the parameters of the object.
const VariablePhiValue & _phi
Shape functions.
Definition: DGKernel.h:98
virtual void precalculateQpResidual(Moose::DGResidualType)
Insertion point for evaluations that depend on qp but are independent of the test functions...
Definition: DGKernel.h:74
virtual Real computeQpResidual(Moose::DGResidualType type)=0
This is the virtual that derived classes should override for computing the residual on neighboring el...
const VariableValue & _u
Holds the current solution at the current quadrature point on the face.
Definition: DGKernel.h:94
OutputTools< Real >::VariableTestGradient VariableTestGradient
Definition: MooseTypes.h:313
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar) override
Computes the element-element off-diagonal Jacobian.
Definition: DGKernel.C:189