https://mooseframework.inl.gov
HFEMTrialJump.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"
13 
27 class HFEMTrialJump : public DGKernel
28 {
29 public:
31 
33 
34 protected:
36  virtual Real computeQpJacobian(Moose::DGJacobianType) override { return 0; }
37  virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar) override;
38 
41 
44 
47 
50 
53 
56 
58  const unsigned int _interior_id;
59 };
virtual Real computeQpJacobian(Moose::DGJacobianType) override
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
Definition: HFEMTrialJump.h:36
const unsigned int _interior_id
Variable id for interior variable.
Definition: HFEMTrialJump.h:58
virtual Real computeQpOffDiagJacobian(Moose::DGJacobianType type, unsigned int jvar) override
This is the virtual that derived classes should override for computing the off-diag Jacobian...
Definition: HFEMTrialJump.C:57
const VariablePhiValue & _test_interior
interior test function at the current side quadrature point
Definition: HFEMTrialJump.h:55
const VariableValue & _interior_neighbor
neighbor's interior solution at the current side quadrature point
Definition: HFEMTrialJump.h:46
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
DGResidualType
Definition: MooseTypes.h:743
virtual Real computeQpResidual(Moose::DGResidualType type) override
This is the virtual that derived classes should override for computing the residual on neighboring el...
Definition: HFEMTrialJump.C:36
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:18
const VariablePhiValue & _phi_interior
interior shape function at the current side quadrature point
Definition: HFEMTrialJump.h:49
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const MooseVariable & _interior_var
Variable whose jump to test against variations in lambda.
Definition: HFEMTrialJump.h:40
OutputTools< Real >::VariableValue VariableValue
Definition: MooseTypes.h:314
DGJacobianType
Definition: MooseTypes.h:749
forward declarations
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const VariableValue & _interior
Current interior solution at the current side quadrature point.
Definition: HFEMTrialJump.h:43
const InputParameters & parameters() const
Get the parameters of the object.
HFEMTrialJump(const InputParameters &parameters)
Definition: HFEMTrialJump.C:23
const VariablePhiValue & _phi_interior_neighbor
neighbor&#39;s shape function at the current side quadrature point
Definition: HFEMTrialJump.h:52
The HFEMTrialJump class computes an L2 product (lambda*, [u])_Gamma, where Gamma is the set of intern...
Definition: HFEMTrialJump.h:27
static InputParameters validParams()
Definition: HFEMTrialJump.C:15