https://mooseframework.inl.gov
HFEMTrialJump.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 "HFEMTrialJump.h"
11 
13 
16 {
18  params.addRequiredCoupledVar("interior_variable", "interior variable to find jumps in");
19  params.addClassDescription("Imposes constraints for HFEM with side-discontinuous variables.");
20  return params;
21 }
22 
24  : DGKernel(parameters),
25  _interior_var(*getVar("interior_variable", 0)),
26  _interior(_is_implicit ? _interior_var.sln() : _interior_var.slnOld()),
27  _interior_neighbor(_is_implicit ? _interior_var.slnNeighbor() : _interior_var.slnOldNeighbor()),
28  _phi_interior(_interior_var.phiFace()),
29  _phi_interior_neighbor(_assembly.phiFaceNeighbor(_interior_var)),
30  _test_interior(_interior_var.phiFace()),
31  _interior_id(coupled("interior_variable"))
32 {
33 }
34 
35 Real
37 {
38  // Use normal vector at qp 0 to make solution depend on geometry
39  // (well, geometry and quadrature rule, with curved boundaries...)
40  // rather than element numbering
41  const Real sign = (_normals[0] > Point()) ? 1 : -1;
42 
43  switch (type)
44  {
45  case Moose::Element:
46  return sign * _test[_i][_qp] * (_interior_neighbor[_qp] - _interior[_qp]);
47  break;
48 
49  case Moose::Neighbor:
50  return 0;
51  break;
52  }
53  return 0;
54 }
55 
56 Real
58 {
59  if (jvar != _interior_id)
60  return 0;
61 
62  const Real sign = (_normals[0] > Point()) ? 1 : -1;
63 
64  switch (type)
65  {
67  return -sign * _test[_i][_qp] * _phi_interior[_j][_qp];
68 
70  return 0;
71 
73  return sign * _test[_i][_qp] * _phi_interior_neighbor[_j][_qp];
74 
76  return 0;
77 
78  default:
79  break;
80  }
81 
82  return 0;
83 }
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 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...
unsigned int _i
Definition: DGKernelBase.h:136
DGResidualType
Definition: MooseTypes.h:744
registerMooseObject("MooseApp", HFEMTrialJump)
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernel.C:27
unsigned int _qp
Definition: DGKernelBase.h:134
unsigned int _j
Definition: DGKernelBase.h:136
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
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
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
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...
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