Loading [MathJax]/extensions/tex2jax.js
www.mooseframework.org
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
DGDiffusion.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "DGDiffusion.h"
11 
12 // MOOSE includes
13 #include "MooseVariableFE.h"
14 
15 #include "libmesh/utility.h"
16 
17 registerMooseObject("MooseApp", DGDiffusion);
18 
20 
23 {
25  // See header file for sigma and epsilon
26  params.addRequiredParam<Real>("sigma", "sigma");
27  params.addRequiredParam<Real>("epsilon", "epsilon");
28  params.addParam<MaterialPropertyName>(
29  "diff", 1., "The diffusion (or thermal conductivity or viscosity) coefficient.");
30  return params;
31 }
32 
34  : DGKernel(parameters),
35  _epsilon(getParam<Real>("epsilon")),
36  _sigma(getParam<Real>("sigma")),
37  _diff(getMaterialProperty<Real>("diff")),
38  _diff_neighbor(getNeighborMaterialProperty<Real>("diff"))
39 {
40 }
41 
42 Real
44 {
45  Real r = 0;
46 
47  const unsigned int elem_b_order = _var.order();
48  const double h_elem =
49  _current_elem->volume() / _current_side_elem->volume() * 1. / Utility::pow<2>(elem_b_order);
50 
51  switch (type)
52  {
53  case Moose::Element:
54  r -= 0.5 * (_diff[_qp] * _grad_u[_qp] * _normals[_qp] +
56  _test[_i][_qp];
57  r += _epsilon * 0.5 * (_u[_qp] - _u_neighbor[_qp]) * _diff[_qp] * _grad_test[_i][_qp] *
58  _normals[_qp];
59  r += _sigma / h_elem * (_u[_qp] - _u_neighbor[_qp]) * _test[_i][_qp];
60  break;
61 
62  case Moose::Neighbor:
63  r += 0.5 * (_diff[_qp] * _grad_u[_qp] * _normals[_qp] +
66  r += _epsilon * 0.5 * (_u[_qp] - _u_neighbor[_qp]) * _diff_neighbor[_qp] *
68  r -= _sigma / h_elem * (_u[_qp] - _u_neighbor[_qp]) * _test_neighbor[_i][_qp];
69  break;
70  }
71 
72  return r;
73 }
74 
75 Real
77 {
78  Real r = 0;
79 
80  const unsigned int elem_b_order = _var.order();
81  const double h_elem =
82  _current_elem->volume() / _current_side_elem->volume() * 1. / Utility::pow<2>(elem_b_order);
83 
84  switch (type)
85  {
87  r -= 0.5 * _diff[_qp] * _grad_phi[_j][_qp] * _normals[_qp] * _test[_i][_qp];
88  r += _epsilon * 0.5 * _phi[_j][_qp] * _diff[_qp] * _grad_test[_i][_qp] * _normals[_qp];
89  r += _sigma / h_elem * _phi[_j][_qp] * _test[_i][_qp];
90  break;
91 
94  r += _epsilon * 0.5 * -_phi_neighbor[_j][_qp] * _diff[_qp] * _grad_test[_i][_qp] *
95  _normals[_qp];
96  r += _sigma / h_elem * -_phi_neighbor[_j][_qp] * _test[_i][_qp];
97  break;
98 
100  r += 0.5 * _diff[_qp] * _grad_phi[_j][_qp] * _normals[_qp] * _test_neighbor[_i][_qp];
101  r += _epsilon * 0.5 * _phi[_j][_qp] * _diff_neighbor[_qp] * _grad_test_neighbor[_i][_qp] *
102  _normals[_qp];
103  r -= _sigma / h_elem * _phi[_j][_qp] * _test_neighbor[_i][_qp];
104  break;
105 
107  r += 0.5 * _diff_neighbor[_qp] * _grad_phi_neighbor[_j][_qp] * _normals[_qp] *
109  r += _epsilon * 0.5 * -_phi_neighbor[_j][_qp] * _diff_neighbor[_qp] *
111  r -= _sigma / h_elem * -_phi_neighbor[_j][_qp] * _test_neighbor[_i][_qp];
112  break;
113  }
114 
115  return r;
116 }
DGKernel::_grad_test_neighbor
const VariableTestGradient & _grad_test_neighbor
Gradient of side shape function.
Definition: DGKernel.h:97
type
MatType type
Definition: PetscDMMoose.C:1477
registerMooseObject
registerMooseObject("MooseApp", DGDiffusion)
Moose::NeighborNeighbor
Definition: MooseTypes.h:646
Moose::DGJacobianType
DGJacobianType
Definition: MooseTypes.h:641
DGKernel::validParams
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernel.C:29
DGDiffusion::_epsilon
Real _epsilon
Definition: DGDiffusion.h:41
MooseObject::type
const std::string & type() const
Get the type of this object.
Definition: MooseObject.h:63
DGKernelBase::_j
unsigned int _j
Definition: DGKernelBase.h:170
Moose::NeighborElement
Definition: MooseTypes.h:645
DGKernelBase::_normals
const MooseArray< Point > & _normals
Normal vectors at the quadrature points.
Definition: DGKernelBase.h:147
DGDiffusion::_sigma
Real _sigma
Definition: DGDiffusion.h:42
InputParameters::addParam
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object.
Definition: InputParameters.h:1198
DGKernel::_grad_phi_neighbor
const VariablePhiGradient & _grad_phi_neighbor
Gradient of side shape function.
Definition: DGKernel.h:93
DGKernel::_u
const VariableValue & _u
Holds the current solution at the current quadrature point on the face.
Definition: DGKernel.h:79
DGKernelBase::_i
unsigned int _i
Definition: DGKernelBase.h:170
DGDiffusion
DG kernel for diffusion.
Definition: DGDiffusion.h:30
defineLegacyParams
defineLegacyParams(DGDiffusion)
DGKernelBase::_current_side_elem
const Elem *& _current_side_elem
Current side element.
Definition: DGKernelBase.h:132
DGKernel::_var
MooseVariable & _var
Variable this kernel operates on.
Definition: DGKernel.h:77
DGDiffusion.h
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
DGKernel::_grad_test
const VariableTestGradient & _grad_test
Gradient of side shape function.
Definition: DGKernel.h:89
DGKernelBase::_current_elem
const Elem *const & _current_elem
Definition: DGKernelBase.h:121
DGKernel::_test_neighbor
const VariableTestValue & _test_neighbor
Side test function.
Definition: DGKernel.h:95
MooseVariableBase::order
Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int.
Definition: MooseVariableBase.C:123
DGKernel::_phi_neighbor
const VariablePhiValue & _phi_neighbor
Side shape function.
Definition: DGKernel.h:91
DGDiffusion::_diff_neighbor
const MaterialProperty< Real > & _diff_neighbor
Definition: DGDiffusion.h:44
Moose::DGResidualType
DGResidualType
Definition: MooseTypes.h:635
Moose::Neighbor
Definition: MooseTypes.h:638
DGDiffusion::DGDiffusion
DGDiffusion(const InputParameters &parameters)
Definition: DGDiffusion.C:33
MooseVariableFE.h
DGKernel::_grad_phi
const VariablePhiGradient & _grad_phi
Gradient of shape function.
Definition: DGKernel.h:85
DGDiffusion::computeQpResidual
virtual Real computeQpResidual(Moose::DGResidualType type) override
This is the virtual that derived classes should override for computing the residual on neighboring el...
Definition: DGDiffusion.C:43
Moose::ElementNeighbor
Definition: MooseTypes.h:644
DGDiffusion::_diff
const MaterialProperty< Real > & _diff
Definition: DGDiffusion.h:43
DGKernel::_u_neighbor
const VariableValue & _u_neighbor
Holds the current solution at the current quadrature point.
Definition: DGKernel.h:99
DGKernelBase::_qp
unsigned int _qp
Definition: DGKernelBase.h:168
DGKernel::_phi
const VariablePhiValue & _phi
Shape functions.
Definition: DGKernel.h:83
Moose::Element
Definition: MooseTypes.h:637
DGKernel::_grad_u_neighbor
const VariableGradient & _grad_u_neighbor
Holds the current solution gradient at the current quadrature point.
Definition: DGKernel.h:101
DGDiffusion::computeQpJacobian
virtual Real computeQpJacobian(Moose::DGJacobianType type) override
This is the virtual that derived classes should override for computing the Jacobian on neighboring el...
Definition: DGDiffusion.C:76
DGKernel::_grad_u
const VariableGradient & _grad_u
Holds the current solution gradient at the current quadrature point on the face.
Definition: DGKernel.h:81
DGKernel::_test
const VariableTestValue & _test
test functions
Definition: DGKernel.h:87
InputParameters::addRequiredParam
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
Definition: InputParameters.h:1176
Moose::ElementElement
Definition: MooseTypes.h:643
DGKernel
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:23
DGDiffusion::validParams
static InputParameters validParams()
Definition: DGDiffusion.C:22