https://mooseframework.inl.gov
InterfaceDiffusion.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 "InterfaceDiffusion.h"
11 
13 
16 {
18  params.addParam<MaterialPropertyName>("D", "D", "The diffusion coefficient.");
19  params.addParam<MaterialPropertyName>(
20  "D_neighbor", "D_neighbor", "The neighboring diffusion coefficient.");
21  params.addClassDescription(
22  "The kernel is utilized to establish flux equivalence on an interface for variables.");
23  return params;
24 }
25 
27  : InterfaceKernel(parameters),
28  _D(getMaterialProperty<Real>("D")),
29  _D_neighbor(getNeighborMaterialProperty<Real>("D_neighbor"))
30 {
31 }
32 
33 Real
35 {
36  Real r = 0;
37 
38  switch (type)
39  {
40  case Moose::Element:
42  break;
43 
44  case Moose::Neighbor:
45  r = _test_neighbor[_i][_qp] * _D[_qp] * _grad_u[_qp] * _normals[_qp];
46  break;
47  }
48 
49  return r;
50 }
51 
52 Real
54 {
55  Real jac = 0;
56 
57  switch (type)
58  {
61  break;
62 
64  jac = _test_neighbor[_i][_qp] * _D[_qp] * _grad_phi[_j][_qp] * _normals[_qp];
65  break;
66 
69  break;
70  }
71 
72  return jac;
73 }
virtual Real computeQpJacobian(Moose::DGJacobianType type) override
Compute jacobians at quadrature points.
InterfaceKernel and VectorInterfaceKernel is responsible for interfacing physics across subdomains...
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
DGResidualType
Definition: MooseTypes.h:743
unsigned int _qp
Current quadrature point.
registerMooseObject("MooseApp", InterfaceDiffusion)
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const MaterialProperty< Real > & _D_neighbor
InterfaceDiffusion(const InputParameters &parameters)
unsigned int _i
Index for test and trial functions.
const TemplateVariableTestValue & _test
Side shape function.
const TemplateVariablePhiGradient & _grad_phi_neighbor
Gradient of side neighbor shape function.
DGJacobianType
Definition: MooseTypes.h:749
const TemplateVariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DG kernel for interfacing diffusion between two variables on adjacent blocks.
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
const TemplateVariablePhiGradient & _grad_phi
Shape function gradient.
virtual Real computeQpResidual(Moose::DGResidualType type) override
Compute residuals at quadrature points.
const TemplateVariableTestValue & _test_neighbor
Side neighbor test function.
const MaterialProperty< Real > & _D
static InputParameters validParams()
const TemplateVariableGradient & _grad_neighbor_value
Coupled neighbor variable gradient.