Line data Source code
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 : 12 : #include "MFEMVectorFEInnerProductIntegralPostprocessor.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMVectorFEInnerProductIntegralPostprocessor); 16 : 17 : InputParameters 18 9170 : MFEMVectorFEInnerProductIntegralPostprocessor::validParams() 19 : { 20 9170 : InputParameters params = MFEMPostprocessor::validParams(); 21 9170 : params += MFEMBlockRestrictable::validParams(); 22 18340 : params.addClassDescription("Calculates the total flux of a vector field through an interface"); 23 36680 : params.addParam<MFEMScalarCoefficientName>( 24 : "coefficient", "1.", "Name of optional scalar coefficient to scale integrand by."); 25 36680 : params.addRequiredParam<VariableName>("primal_variable", 26 : "Name of the first vector variable in the inner product."); 27 27510 : params.addRequiredParam<VariableName>("dual_variable", 28 : "Name of the second vector variable in the inner product."); 29 9170 : return params; 30 0 : } 31 : 32 30 : MFEMVectorFEInnerProductIntegralPostprocessor::MFEMVectorFEInnerProductIntegralPostprocessor( 33 30 : const InputParameters & parameters) 34 : : MFEMPostprocessor(parameters), 35 30 : MFEMBlockRestrictable(parameters, getMFEMProblem().mesh().getMFEMParMesh()), 36 30 : _primal_var(getMFEMProblem().getProblemData().gridfunctions.GetRef( 37 60 : getParam<VariableName>("primal_variable"))), 38 30 : _dual_var(getMFEMProblem().getProblemData().gridfunctions.GetRef( 39 60 : getParam<VariableName>("dual_variable"))), 40 60 : _scalar_coef(getScalarCoefficient("coefficient")), 41 30 : _dual_var_coef(&_dual_var), 42 30 : _scaled_dual_var_coef(_scalar_coef, _dual_var_coef), 43 90 : _subdomain_integrator(_primal_var.ParFESpace()) 44 : { 45 30 : if (isSubdomainRestricted()) 46 : { 47 60 : _subdomain_integrator.AddDomainIntegrator( 48 30 : new mfem::VectorFEDomainLFIntegrator(_scaled_dual_var_coef), getSubdomainMarkers()); 49 : } 50 : else 51 : { 52 0 : _subdomain_integrator.AddDomainIntegrator( 53 0 : new mfem::VectorFEDomainLFIntegrator(_scaled_dual_var_coef)); 54 : } 55 30 : } 56 : 57 : void 58 51 : MFEMVectorFEInnerProductIntegralPostprocessor::execute() 59 : { 60 51 : _subdomain_integrator.Assemble(); 61 51 : _integral = _subdomain_integrator(_primal_var); 62 51 : } 63 : 64 : PostprocessorValue 65 51 : MFEMVectorFEInnerProductIntegralPostprocessor::getValue() const 66 : { 67 51 : return _integral; 68 : } 69 : 70 : #endif