https://mooseframework.inl.gov
MFEMVectorFEInnerProductIntegralPostprocessor.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 #include "MFEMProblem.h"
14 
16 
19 {
22  params.addClassDescription("Calculates the total flux of a vector field through an interface");
23  params.addParam<MFEMScalarCoefficientName>(
24  "coefficient", "1.", "Name of optional scalar coefficient to scale integrand by.");
25  params.addRequiredParam<VariableName>("primal_variable",
26  "Name of the first vector variable in the inner product.");
27  params.addRequiredParam<VariableName>("dual_variable",
28  "Name of the second vector variable in the inner product.");
29  return params;
30 }
31 
33  const InputParameters & parameters)
34  : MFEMPostprocessor(parameters),
35  MFEMBlockRestrictable(parameters, getMFEMProblem().mesh().getMFEMParMesh()),
36  _primal_var(getMFEMProblem().getProblemData().gridfunctions.GetRef(
37  getParam<VariableName>("primal_variable"))),
38  _dual_var(getMFEMProblem().getProblemData().gridfunctions.GetRef(
39  getParam<VariableName>("dual_variable"))),
40  _scalar_coef(getScalarCoefficient("coefficient")),
41  _dual_var_coef(&_dual_var),
42  _scaled_dual_var_coef(_scalar_coef, _dual_var_coef),
43  _subdomain_integrator(_primal_var.ParFESpace())
44 {
46  {
47  _subdomain_integrator.AddDomainIntegrator(
48  new mfem::VectorFEDomainLFIntegrator(_scaled_dual_var_coef), getSubdomainMarkers());
49  }
50  else
51  {
52  _subdomain_integrator.AddDomainIntegrator(
53  new mfem::VectorFEDomainLFIntegrator(_scaled_dual_var_coef));
54  }
55 }
56 
57 void
59 {
60  _subdomain_integrator.Assemble();
62 }
63 
66 {
67  return _integral;
68 }
69 
70 #endif
mfem::Array< int > & getSubdomainMarkers()
MeshBase & mesh
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
bool isSubdomainRestricted()
Returns a bool indicating if the object is restricted to a subset of subdomains.
Postprocessor for MFEM results.
Real PostprocessorValue
various MOOSE typedefs
Definition: MooseTypes.h:202
virtual PostprocessorValue getValue() const override final
Return the last evaluated integral value.
registerMooseObject("MooseApp", MFEMVectorFEInnerProductIntegralPostprocessor)
static InputParameters validParams()
Compute the integral of the innter product between two MFEM vector FE variables, scaled by an optiona...
static InputParameters validParams()
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...
Base class for construction of an object that is restricted to a subset of subdomains of the problem ...