https://mooseframework.inl.gov
MassMatrixDGKernel.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 "MassMatrixDGKernel.h"
11 
12 registerMooseObject("NavierStokesApp", MassMatrixDGKernel);
13 
16 {
18  params.addClassDescription("Computes a finite element mass matrix meant for use in "
19  "preconditioning schemes which require one");
20  params.addParam<Real>("density", 1, "Optional density for scaling the computed mass.");
21  params.set<MultiMooseEnum>("vector_tags") = "";
22  params.set<MultiMooseEnum>("matrix_tags") = "";
23  params.suppressParameter<MultiMooseEnum>("vector_tags");
24  params.suppressParameter<std::vector<TagName>>("extra_vector_tags");
25  params.suppressParameter<std::vector<TagName>>("absolute_value_vector_tags");
26  params.set<bool>("matrix_only") = true;
27  return params;
28 }
29 
31  : DGKernel(parameters), _density(getParam<Real>("density"))
32 {
33  if (!isParamValid("matrix_tags") && !isParamValid("extra_matrix_tags"))
34  mooseError("One of 'matrix_tags' or 'extra_matrix_tags' must be provided");
35 }
36 
37 Real
39 {
40  mooseAssert(false, "should never be called");
41  return 0;
42 }
43 
44 Real
46 {
47  Real jac = 0;
48 
49  switch (type)
50  {
52  jac = _test[_i][_qp] * _density * _phi[_j][_qp];
53  break;
54 
55  default:
56  jac = 0;
57  break;
58  }
59 
60  return jac;
61 }
This class can be used to compute a mass matrix for facet unknowns, e.g.
MassMatrixDGKernel(const InputParameters &parameters)
static InputParameters validParams()
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
T & set(const std::string &name, bool quiet_mode=false)
unsigned int _i
DGResidualType
ElementElement
void suppressParameter(const std::string &name)
bool isParamValid(const std::string &name) const
static InputParameters validParams()
unsigned int _qp
unsigned int _j
const std::string & type() const
const VariableTestValue & _test
DGJacobianType
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("NavierStokesApp", MassMatrixDGKernel)
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
const VariablePhiValue & _phi
virtual Real computeQpJacobian(const Moose::DGJacobianType type) override
virtual Real computeQpResidual(Moose::DGResidualType) override