https://mooseframework.inl.gov
MFEMVectorFEMassKernel.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 #ifdef MFEM_ENABLED
11 
12 #include "MFEMVectorFEMassKernel.h"
13 #include "MFEMProblem.h"
14 
16 
19 {
21  params.addClassDescription("Adds the domain integrator to an MFEM problem for the bilinear form "
22  "$(k \\vec u, \\vec v)_\\Omega$ "
23  "arising from the weak form of the mass operator "
24  "$k \\vec u$.");
25  params.addParam<MFEMScalarCoefficientName>(
26  "coefficient", "1.", "Name of property k to multiply the integrator by");
27  return params;
28 }
29 
31  : MFEMKernel(parameters),
32  // FIXME: The MFEM bilinear form can also handle vector and matrix
33  // coefficients, so ideally we'd handle all three too.
34  _coef(getScalarCoefficient(getParam<MFEMScalarCoefficientName>("coefficient")))
35 {
36 }
37 
38 mfem::BilinearFormIntegrator *
40 {
41  return new mfem::VectorFEMassIntegrator(_coef);
42 }
43 
44 #endif
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual mfem::BilinearFormIntegrator * createBFIntegrator() override
MFEMVectorFEMassKernel(const InputParameters &parameters)
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...
static InputParameters validParams()
static InputParameters validParams()
Definition: MFEMKernel.C:19
registerMooseObject("MooseApp", MFEMVectorFEMassKernel)