https://mooseframework.inl.gov
HFEMDiffusion.C
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 #include "HFEMDiffusion.h"
11 
13 
16 {
18  params.addClassDescription("Imposes the constraints on internal sides with HFEM.");
19  return params;
20 }
21 
22 HFEMDiffusion::HFEMDiffusion(const InputParameters & parameters) : DGLowerDKernel(parameters) {}
23 
24 Real
26 {
27  switch (type)
28  {
29  case Moose::Element:
30  return -_lambda[_qp] * _test[_i][_qp];
31 
32  case Moose::Neighbor:
33  return _lambda[_qp] * _test_neighbor[_i][_qp];
34  }
35  return 0;
36 }
37 
38 Real
40 {
41  return (_u_neighbor[_qp] - _u[_qp]) * _test_lambda[_i][_qp];
42 }
43 
45 
46 Real
48 {
49  switch (type)
50  {
52  return -_test_lambda[_i][_qp] * _phi[_j][_qp];
53 
55  return _test_lambda[_i][_qp] * _phi_neighbor[_j][_qp];
56 
58  return -_phi_lambda[_j][_qp] * _test[_i][_qp];
59 
61  return _phi_lambda[_j][_qp] * _test_neighbor[_i][_qp];
62 
63  default:
64  break;
65  }
66 
67  return 0;
68 }
const VariableValue & _u_neighbor
Holds the current solution at the current quadrature point.
Definition: DGKernel.h:114
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...
unsigned int _i
Definition: DGKernelBase.h:136
DGResidualType
Definition: MooseTypes.h:743
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
const VariableTestValue & _test_neighbor
Side test function.
Definition: DGKernel.h:110
unsigned int _qp
Definition: DGKernelBase.h:134
const VariableValue & _lambda
Holds the current solution at the current quadrature point on the face.
const VariablePhiValue & _phi_lambda
Shape functions.
unsigned int _j
Definition: DGKernelBase.h:136
const VariableTestValue & _test_lambda
test functions
static InputParameters validParams()
Definition: HFEMDiffusion.C:15
virtual Real computeQpJacobian(Moose::DGJacobianType type) override
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
Definition: HFEMDiffusion.C:44
registerMooseObject("MooseApp", HFEMDiffusion)
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
HFEMDiffusion(const InputParameters &parameters)
Definition: HFEMDiffusion.C:22
virtual Real computeQpResidual(Moose::DGResidualType type) override
This is the virtual that derived classes should override for computing the residual on neighboring el...
Definition: HFEMDiffusion.C:25
const VariableTestValue & _test
test functions
Definition: DGKernel.h:102
DGJacobianType
Definition: MooseTypes.h:749
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType type) override
Computes one of the five pieces of Jacobian involving lower-d at quadrature points.
Definition: HFEMDiffusion.C:47
ConstraintJacobianType
Definition: MooseTypes.h:796
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
const VariablePhiValue & _phi
Shape functions.
Definition: DGKernel.h:98
const VariableValue & _u
Holds the current solution at the current quadrature point on the face.
Definition: DGKernel.h:94
virtual Real computeLowerDQpResidual() override
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
Definition: HFEMDiffusion.C:39