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 MFEM_ENABLED 11 : 12 : #include "MFEMDivAux.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMDivAux); 16 : 17 : /* 18 : Class to set an L2 auxvariable to be the divergence of a H(div) vector variable. 19 : */ 20 : InputParameters 21 8646 : MFEMDivAux::validParams() 22 : { 23 8646 : InputParameters params = MFEMAuxKernel::validParams(); 24 8646 : params.addClassDescription( 25 : "Calculates the divergence of an H(div) conforming RT source variable and stores the result" 26 : " on an L2 conforming result auxvariable"); 27 8646 : params.addRequiredParam<VariableName>("source", 28 : "Vector H(div) MFEMVariable to take the divergence of."); 29 8646 : params.addParam<mfem::real_t>("scale_factor", 1.0, "Factor to scale result auxvariable by."); 30 8646 : return params; 31 0 : } 32 : 33 8 : MFEMDivAux::MFEMDivAux(const InputParameters & parameters) 34 : : MFEMAuxKernel(parameters), 35 8 : _source_var_name(getParam<VariableName>("source")), 36 8 : _source_var(*getMFEMProblem().getProblemData().gridfunctions.Get(_source_var_name)), 37 8 : _scale_factor(getParam<mfem::real_t>("scale_factor")), 38 16 : _div(_source_var.ParFESpace(), _result_var.ParFESpace()) 39 : { 40 8 : _div.Assemble(); 41 8 : _div.Finalize(); 42 8 : } 43 : 44 : // Computes the auxvariable. 45 : void 46 8 : MFEMDivAux::execute() 47 : { 48 8 : _result_var = 0.0; 49 8 : _div.AddMult(_source_var, _result_var, _scale_factor); 50 8 : } 51 : 52 : #endif