https://mooseframework.inl.gov
FVMassMatrix.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 "FVMassMatrix.h"
11 #include "Function.h"
12 
14 
17 {
19  params.addClassDescription(
20  "Computes a 'mass matrix', which will just be a diagonal matrix for the finite volume "
21  "method, meant for use in preconditioning schemes which require one");
22  params.addParam<MooseFunctorName>("density", 1.0, "Optional density weighting functor");
23  params.set<MultiMooseEnum>("vector_tags") = "";
24  params.set<MultiMooseEnum>("matrix_tags") = "";
25  params.suppressParameter<MultiMooseEnum>("vector_tags");
26  params.suppressParameter<std::vector<TagName>>("extra_vector_tags");
27  params.suppressParameter<std::vector<TagName>>("absolute_value_vector_tags");
28  params.set<bool>("matrix_only") = true;
29  return params;
30 }
31 
33  : FVElementalKernel(parameters), _density(getFunctor<Real>("density"))
34 {
35  if (!isParamValid("matrix_tags") && !isParamValid("extra_matrix_tags"))
36  mooseError("One of 'matrix_tags' or 'extra_matrix_tags' must be provided");
37 }
38 
39 void
41 {
42 }
43 
44 ADReal
46 {
47  const auto elem = makeElemArg(_current_elem);
48  const auto state = determineState();
49  return _density(elem, state) * _var(elem, state);
50 }
ADReal computeQpResidual() override
This is the primary function that must be implemented for flux kernel terms.
Definition: FVMassMatrix.C:45
Moose::StateArg determineState() const
Create a functor state argument that corresponds to the implicit state of this object.
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
FVMassMatrix(const InputParameters &parameters)
Definition: FVMassMatrix.C:32
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
void suppressParameter(const std::string &name)
This method suppresses an inherited parameter so that it isn&#39;t required or valid in the derived class...
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
Helper method to create an elemental argument for a functor that includes whether to perform skewness...
FVElemental is used for calculating residual contributions from volume integral terms of a PDE where ...
const Elem *const & _current_elem
static InputParameters validParams()
registerMooseObject("MooseApp", FVMassMatrix)
const Moose::Functor< Real > & _density
Elemental weighting functor.
Definition: FVMassMatrix.h:31
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Computes a &#39;mass matrix&#39;, which will just be a diagonal matrix for the finite volume method...
Definition: FVMassMatrix.h:18
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
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...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
static InputParameters validParams()
Definition: FVMassMatrix.C:16
MooseVariableFV< Real > & _var
virtual void computeResidual() override
Usually you should not override these functions - they have some tricky stuff in them that you don&#39;t ...
Definition: FVMassMatrix.C:40