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 "MFEMVectorBoundaryFluxIntegralPostprocessor.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMVectorBoundaryFluxIntegralPostprocessor); 16 : 17 : InputParameters 18 2124 : MFEMVectorBoundaryFluxIntegralPostprocessor::validParams() 19 : { 20 2124 : InputParameters params = MFEMPostprocessor::validParams(); 21 2124 : params += MFEMBoundaryRestrictable::validParams(); 22 4248 : params.addClassDescription( 23 : "Calculates the integral of the flux of a vector variable across a boundary."); 24 8496 : params.addParam<MFEMScalarCoefficientName>( 25 : "coefficient", "1.", "Name of optional scalar coefficient to scale integrand by."); 26 6372 : params.addRequiredParam<VariableName>("variable", "Name of the vector variable."); 27 2124 : return params; 28 0 : } 29 : 30 13 : MFEMVectorBoundaryFluxIntegralPostprocessor::MFEMVectorBoundaryFluxIntegralPostprocessor( 31 13 : const InputParameters & parameters) 32 : : MFEMPostprocessor(parameters), 33 : MFEMBoundaryRestrictable( 34 39 : parameters, getMFEMProblem().getMFEMVariableMesh(getParam<VariableName>("variable"))), 35 26 : _var(*getMFEMProblem().getGridFunction(getParam<VariableName>("variable"))), 36 26 : _scalar_coef(getScalarCoefficient("coefficient")), 37 13 : _var_coef(&_var), 38 13 : _rt_fec(_var.ParFESpace()->GetMaxElementOrder(), getMesh().Dimension()), 39 13 : _rt_vector_fespace(const_cast<mfem::ParMesh *>(&getMesh()), &_rt_fec), 40 13 : _rt_var(&_rt_vector_fespace), 41 52 : _boundary_integrator(&_rt_vector_fespace) 42 : { 43 13 : _rt_var = 0; 44 26 : _boundary_integrator.AddBoundaryIntegrator( 45 13 : new mfem::VectorFEBoundaryFluxLFIntegrator(_scalar_coef), getBoundaryMarkers()); 46 13 : } 47 : 48 : void 49 13 : MFEMVectorBoundaryFluxIntegralPostprocessor::execute() 50 : { 51 13 : _rt_var.ProjectBdrCoefficientNormal(_var_coef, getBoundaryMarkers()); 52 13 : _boundary_integrator.Assemble(); 53 13 : _integral = _boundary_integrator(_rt_var); 54 13 : } 55 : 56 : PostprocessorValue 57 13 : MFEMVectorBoundaryFluxIntegralPostprocessor::getValue() const 58 : { 59 13 : return _integral; 60 : } 61 : 62 : #endif