https://mooseframework.inl.gov
HFEMTestJump.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 "HFEMTestJump.h"
11 
13 
16 {
18  params.addRequiredCoupledVar("side_variable", "side variable to use as Lagrange multiplier");
19  params.addClassDescription("Imposes constraints for HFEM with side-discontinuous variables.");
20  return params;
21 }
22 
24  : DGKernel(parameters),
25  _lambda_var(*getVar("side_variable", 0)),
26  _lambda(_is_implicit ? _lambda_var.sln() : _lambda_var.slnOld()),
27  _phi_lambda(_lambda_var.phiFace()),
28  _test_lambda(_lambda_var.phiFace()),
29  _lambda_id(coupled("side_variable"))
30 {
31 }
32 
33 Real
35 {
36  // Use normal vector at qp 0 to make solution depend on geometry
37  // (well, geometry and quadrature rule, with curved boundaries...)
38  // rather than element numbering
39  const Real sign = (_normals[0] > Point()) ? 1 : -1;
40 
41  switch (type)
42  {
43  case Moose::Element:
44  return -sign * _lambda[_qp] * _test[_i][_qp];
45  break;
46 
47  case Moose::Neighbor:
48  return sign * _lambda[_qp] * _test_neighbor[_i][_qp];
49  break;
50  }
51  return 0;
52 }
53 
54 Real
56 {
57  if (jvar != _lambda_id)
58  return 0;
59 
60  const Real sign = (_normals[0] > Point()) ? 1 : -1;
61 
62  switch (type)
63  {
65  return -sign * _phi_lambda[_j][_qp] * _test[_i][_qp];
66 
68  return sign * _phi_lambda[_j][_qp] * _test_neighbor[_i][_qp];
69 
71  return 0;
72 
74  return 0;
75 
76  default:
77  break;
78  }
79 
80  return 0;
81 }
registerMooseObject("MooseApp", HFEMTestJump)
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:744
The HFEMTestJump class computes an L2 product (lambda, [u*])_Gamma, where Gamma is the set of interna...
Definition: HFEMTestJump.h:27
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
unsigned int _qp
Definition: DGKernelBase.h:134
unsigned int _j
Definition: DGKernelBase.h:136
const unsigned int _lambda_id
Variable id for lambda.
Definition: HFEMTestJump.h:52
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:18
virtual Real computeQpResidual(Moose::DGResidualType type) override
This is the virtual that derived classes should override for computing the residual on neighboring el...
Definition: HFEMTestJump.C:34
HFEMTestJump(const InputParameters &parameters)
Definition: HFEMTestJump.C:23
T sign(T x)
Definition: MathUtils.h:84
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
Definition: DGKernelBase.h:113
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
const VariableTestValue & _test
test functions
Definition: DGKernel.h:102
DGJacobianType
Definition: MooseTypes.h:750
static InputParameters validParams()
Definition: HFEMTestJump.C:15
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 VariableValue & _lambda
Current lambda solution at the current side quadrature point.
Definition: HFEMTestJump.h:43
const VariablePhiValue & _phi_lambda
lambda shape function at the current side quadrature point
Definition: HFEMTestJump.h:46
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: HFEMTestJump.C:55