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 9406 : MFEMSumAux::validParams() 19 : { 20 9406 : InputParameters params = MFEMAuxKernel::validParams(); 21 18812 : params.addClassDescription( 22 : "Calculates the sum of two variables sharing an FE space, each optionally scaled by a real " 23 : "constant, and stores the result in a third."); 24 37624 : params.addRequiredParam<std::vector<VariableName>>("source_variables", 25 : "The names of MFEM variables to sum over"); 26 28218 : params.addParam<std::vector<mfem::real_t>>( 27 : "scale_factors", "The factors to scale each MFEM variable by during summation"); 28 9406 : return params; 29 0 : } 30 : 31 29 : MFEMSumAux::MFEMSumAux(const InputParameters & parameters) 32 : : MFEMAuxKernel(parameters), 33 58 : _var_names(getParam<std::vector<VariableName>>("source_variables")), 34 122 : _scale_factors(parameters.isParamValid("scale_factors") 35 29 : ? getParam<std::vector<mfem::real_t>>("scale_factors") 36 110 : : std::vector<mfem::real_t>(_var_names.size(), 1.0)) 37 : { 38 29 : 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 81 : for (const auto & var_name : _var_names) 43 : { 44 : const mfem::ParGridFunction * gf = 45 56 : getMFEMProblem().getProblemData().gridfunctions.Get(var_name); 46 56 : if (gf->ParFESpace() == _result_var.ParFESpace()) 47 54 : _summed_vars.push_back(gf); 48 : else 49 10 : paramError("source_variables", 50 : "The MFEM variable ", 51 : var_name, 52 : " being summed has a different FESpace from ", 53 2 : _result_var_name); 54 : } 55 37 : } 56 : 57 : void 58 25 : MFEMSumAux::execute() 59 : { 60 : // result = sum_i (_scale_factor_i * _summed_var_i) 61 25 : _result_var = 0.0; 62 77 : for (const auto i : index_range(_summed_vars)) 63 52 : _result_var.Add(_scale_factors[i], *_summed_vars[i]); 64 25 : } 65 : 66 : #endif