www.mooseframework.org
AEFVKernel.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 "AEFVKernel.h"
11 
13 
16 {
18  params.addClassDescription(
19  "A dgkernel for the advection equation using a cell-centered finite volume method.");
20  MooseEnum component("concentration");
21  params.addParam<MooseEnum>("component", component, "Choose one of the equations");
22  params.addRequiredCoupledVar("u", "Name of the variable to use");
23  params.addRequiredParam<UserObjectName>("flux", "Name of the internal side flux object to use");
24  return params;
25 }
26 
28  : DGKernel(parameters),
29  _component(getParam<MooseEnum>("component")),
30  _uc1(coupledValue("u")),
31  _uc2(coupledNeighborValue("u")),
32  _u1(getMaterialProperty<Real>("u")),
33  _u2(getNeighborMaterialProperty<Real>("u")),
34  _flux(getUserObject<InternalSideFluxBase>("flux"))
35 {
36 }
37 
39 
40 Real
42 {
43  // assemble the input vectors, which are
44  // the reconstructed linear monomial
45  // extrapolated at side center from the current and neighbor elements
46  std::vector<Real> uvec1 = {_u1[_qp]};
47  std::vector<Real> uvec2 = {_u2[_qp]};
48 
49  // calculate the flux
50  const auto & flux = _flux.getFlux(
51  _current_side, _current_elem->id(), _neighbor_elem->id(), uvec1, uvec2, _normals[_qp]);
52 
53  // distribute the contribution to the current and neighbor elements
54  switch (type)
55  {
56  case Moose::Element:
57  return flux[_component] * _test[_i][_qp];
58 
59  case Moose::Neighbor:
60  return -flux[_component] * _test_neighbor[_i][_qp];
61  }
62 
63  return 0.0;
64 }
65 
66 Real
68 {
69  // assemble the input vectors, which are
70  // the constant monomial from the current and neighbor elements
71  std::vector<Real> uvec1 = {_uc1[_qp]};
72  std::vector<Real> uvec2 = {_uc2[_qp]};
73 
74  // calculate the Jacobian matrices
75  const auto & fjac1 = _flux.getJacobian(Moose::Element,
77  _current_elem->id(),
78  _neighbor_elem->id(),
79  uvec1,
80  uvec2,
81  _normals[_qp]);
82 
83  const auto & fjac2 = _flux.getJacobian(Moose::Neighbor,
85  _current_elem->id(),
86  _neighbor_elem->id(),
87  uvec1,
88  uvec2,
89  _normals[_qp]);
90 
91  // distribute the contribution to the current and neighbor elements
92  switch (type)
93  {
95  return fjac1(_component, _component) * _phi[_j][_qp] * _test[_i][_qp];
96 
98  return fjac2(_component, _component) * _phi_neighbor[_j][_qp] * _test[_i][_qp];
99 
101  return -fjac1(_component, _component) * _phi[_j][_qp] * _test_neighbor[_i][_qp];
102 
104  return -fjac2(_component, _component) * _phi_neighbor[_j][_qp] * _test_neighbor[_i][_qp];
105  }
106 
107  return 0.0;
108 }
NeighborElement
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static const std::string component
Definition: NS.h:138
registerMooseObject("RdgApp", AEFVKernel)
const VariablePhiValue & _phi_neighbor
unsigned int _i
DGResidualType
virtual Real computeQpJacobian(Moose::DGJacobianType type)
Definition: AEFVKernel.C:67
const InternalSideFluxBase & _flux
flux user object
Definition: AEFVKernel.h:65
ElementElement
static InputParameters validParams()
Definition: AEFVKernel.C:15
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
const VariableTestValue & _test_neighbor
const MaterialProperty< Real > & _u2
Definition: AEFVKernel.h:62
unsigned int _qp
unsigned int _j
const MaterialProperty< Real > & _u1
extrapolated variable values at side center
Definition: AEFVKernel.h:61
A dgkernel for the advection equation using a cell-centered finite volume method. ...
Definition: AEFVKernel.h:38
ElementNeighbor
virtual const DenseMatrix< Real > & getJacobian(Moose::DGResidualType type, unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave) const
Get the Jacobian matrix.
const std::string & type() const
const MooseArray< Point > & _normals
const VariableValue & _uc1
piecewise constant variable values in cells
Definition: AEFVKernel.h:57
AEFVKernel(const InputParameters &parameters)
Definition: AEFVKernel.C:27
const VariableValue & _uc2
Definition: AEFVKernel.h:58
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
const VariableTestValue & _test
DGJacobianType
MooseEnum _component
choose an equation
Definition: AEFVKernel.h:51
virtual Real computeQpResidual(Moose::DGResidualType type)
Definition: AEFVKernel.C:41
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual const std::vector< Real > & getFlux(unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave) const
Get the flux vector.
A base class for computing and caching internal side flux.
const Elem *const & _current_elem
const unsigned int & _current_side
void addClassDescription(const std::string &doc_string)
const Elem *const & _neighbor_elem
const VariablePhiValue & _phi
NeighborNeighbor
virtual ~AEFVKernel()
Definition: AEFVKernel.C:38