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 "MFEMSumAux.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMSumAux); 16 : 17 : InputParameters 18 2206 : MFEMSumAux::validParams() 19 : { 20 2206 : InputParameters params = MFEMAuxKernel::validParams(); 21 4412 : params.addClassDescription( 22 : "Calculates the sum of an arbitrary number of variables sharing an FE space, each " 23 : "optionally scaled by a real constant, and stores the result in an auxiliary variable."); 24 8824 : MFEMExecutedObject::addRequiredDependencyParam<std::vector<VariableName>>( 25 : params, "source_variables", "The names of the MFEM variables to sum"); 26 6618 : params.addParam<std::vector<mfem::real_t>>( 27 : "scale_factors", "The factors to scale each MFEM variable by during summation"); 28 2206 : return params; 29 0 : } 30 : 31 56 : MFEMSumAux::MFEMSumAux(const InputParameters & parameters) 32 : : MFEMAuxKernel(parameters), 33 112 : _var_names(getParam<std::vector<VariableName>>("source_variables")), 34 230 : _scale_factors(parameters.isParamValid("scale_factors") 35 56 : ? getParam<std::vector<mfem::real_t>>("scale_factors") 36 218 : : std::vector<mfem::real_t>(_var_names.size(), 1.0)) 37 : { 38 56 : if (_var_names.size() != _scale_factors.size()) 39 6 : paramError("scale_factors", 40 : "Number of MFEM variables to sum over is different from the number of provided " 41 : "scale factors."); 42 162 : for (const auto & var_name : _var_names) 43 : { 44 110 : const mfem::ParGridFunction * gf = getMFEMProblem().getGridFunction(var_name).get(); 45 110 : if (gf->ParFESpace() == _result_var.ParFESpace()) 46 108 : _summed_vars.push_back(gf); 47 : else 48 10 : paramError("source_variables", 49 : "The MFEM variable ", 50 : var_name, 51 : " being summed has a different FESpace from ", 52 2 : _result_var_name); 53 : } 54 64 : } 55 : 56 : void 57 52 : MFEMSumAux::execute() 58 : { 59 : // result = sum_i (_scale_factor_i * _summed_var_i) 60 52 : _result_var = 0.0; 61 158 : for (const auto i : index_range(_summed_vars)) 62 106 : _result_var.Add(_scale_factors[i], *_summed_vars[i]); 63 52 : } 64 : 65 : #endif