https://mooseframework.inl.gov
MFEMMixedBilinearFormKernel.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 MOOSE_MFEM_ENABLED
11 
13 
16 {
18  params.addClassDescription(
19  "Base class for mixed bilinear form kernels, allowing different trial and test variables.");
20  params.addParam<VariableName>(
21  "trial_variable",
22  "The trial variable this kernel is acting on and which will be solved for. If empty "
23  "(default), it will be the same as the test variable.");
24  params.addParam<bool>(
25  "transpose", false, "If true, adds the transpose of the integrator to the system instead.");
26  return params;
27 }
28 
30  : MFEMKernel(parameters),
31  _trial_var_name(isParamValid("trial_variable") ? getParam<VariableName>("trial_variable")
32  : _test_var_name),
33  _transpose(getParam<bool>("transpose"))
34 {
35 }
36 
37 const VariableName &
39 {
40  return _trial_var_name;
41 }
42 
43 mfem::BilinearFormIntegrator *
45 {
46  return _transpose ? new mfem::TransposeIntegrator(createMBFIntegrator()) : createMBFIntegrator();
47 }
48 #endif
MFEMMixedBilinearFormKernel(const InputParameters &parameters)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual const VariableName & getTrialVariableName() const override
Get name of the trial variable (gridfunction) the kernel acts on.
const VariableName _trial_var_name
Name of the trial variable that the kernel is applied to.
virtual mfem::BilinearFormIntegrator * createBFIntegrator() override
We override this to optionally transpose the mixed bilinear form integrator.
static InputParameters validParams()
virtual mfem::BilinearFormIntegrator * createMBFIntegrator()
Create MFEM mixed bilinear form integrator. Ownership managed by the caller.
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...
bool _transpose
Bool controlling whether to add the transpose of the integrator to the system.
Class to construct an MFEM integrator to apply to the equation system.
Definition: MFEMKernel.h:21
static InputParameters validParams()
Definition: MFEMKernel.C:19