Skip to content
10 changes: 5 additions & 5 deletions framework/doc/content/source/mfem/auxkernels/MFEMSumAux.md
Original file line number Diff line number Diff line change
Expand Up @@ -8,14 +8,14 @@

## Overview

AuxKernel for calculating the sum of two MFEM variables that are defined on the same finite element
space, and storing the result in a third. Both variables may be (optionally) scaled by a real scalar
constant prior to addition.
AuxKernel for calculating the sum of two or more MFEM variables that are defined on the same finite
element space, and storing the result in another. All variables may be (optionally) scaled by a real
scalar constant prior to addition.

!equation
u = \lambda_1 v_1 + \lambda_2 v_2.
u = \sum_i \lambda_i v_i.

where $u, v_1, v_2$ are defined on the same FE space, and $\lambda_1, \lambda_2$ are real scalar
where $u$ and $\{v_i\}$ are defined on the same FE space, and $\{\lambda_i\}$ are real scalar
constants.

## Example Input File Syntax
Expand Down
11 changes: 4 additions & 7 deletions framework/include/mfem/auxkernels/MFEMSumAux.h
Original file line number Diff line number Diff line change
Expand Up @@ -32,14 +32,11 @@ class MFEMSumAux : public MFEMAuxKernel

protected:
// Names of input MFEMVariables to sum.
const VariableName & _v1_var_name;
const VariableName & _v2_var_name;
// Reference to input variable gridfunctions.
const mfem::ParGridFunction & _v1_var;
const mfem::ParGridFunction & _v2_var;
const std::vector<VariableName> & _var_names;
// Scalar factors to multiply the input variables by.
const mfem::real_t _lambda1;
const mfem::real_t _lambda2;
const std::vector<mfem::real_t> _scale_factors;
// References to input variable gridfunctions.
std::vector<const mfem::ParGridFunction *> _summed_vars;
};

#endif
45 changes: 31 additions & 14 deletions framework/src/mfem/auxkernels/MFEMSumAux.C
Original file line number Diff line number Diff line change
Expand Up @@ -21,31 +21,48 @@ MFEMSumAux::validParams()
params.addClassDescription(
"Calculates the sum of two variables sharing an FE space, each optionally scaled by a real "
"constant, and stores the result in a third.");
params.addRequiredParam<VariableName>("first_source_variable", "First variable to sum.");
params.addRequiredParam<VariableName>("second_source_variable", "Second variable to sum.");
params.addParam<mfem::real_t>(
"first_scale_factor", 1.0, "Factor to scale the first variable by prior to sum.");
params.addParam<mfem::real_t>(
"second_scale_factor", 1.0, "Factor to scale the second variable by prior to sum.");
params.addRequiredParam<std::vector<VariableName>>("source_variables",
"The names of MFEM variables to sum over");
params.addParam<std::vector<mfem::real_t>>(
"scale_factors", "The factors to scale each MFEM variable by during summation");
return params;
}

MFEMSumAux::MFEMSumAux(const InputParameters & parameters)
: MFEMAuxKernel(parameters),
_v1_var_name(getParam<VariableName>("first_source_variable")),
_v2_var_name(getParam<VariableName>("second_source_variable")),
_v1_var(*getMFEMProblem().getProblemData().gridfunctions.Get(_v1_var_name)),
_v2_var(*getMFEMProblem().getProblemData().gridfunctions.Get(_v2_var_name)),
_lambda1(getParam<mfem::real_t>("first_scale_factor")),
_lambda2(getParam<mfem::real_t>("second_scale_factor"))
_var_names(getParam<std::vector<VariableName>>("source_variables")),
_scale_factors(parameters.isParamValid("scale_factors")
? getParam<std::vector<mfem::real_t>>("scale_factors")
: std::vector<mfem::real_t>(_var_names.size(), 1.0))
{
if (_var_names.size() != _scale_factors.size())
paramError("scale_factors",
"Number of MFEM variables to sum over is different from the number of provided "
"scale factors.");
for (const auto & var_name : _var_names)
{
const mfem::ParGridFunction * gf =
getMFEMProblem().getProblemData().gridfunctions.Get(var_name);
if (gf->ParFESpace() == _result_var.ParFESpace())
_summed_vars.push_back(gf);
else
paramError("source_variables",
"The MFEM variable ",
var_name,
" being summed has a different FESpace from ",
_result_var_name);
}
}

void
MFEMSumAux::execute()
{
// result = lambda1 * v1 + lambda2 * v2
add(_lambda1, _v1_var, _lambda2, _v2_var, _result_var);
// result = sum_i (_scale_factor_i * _summed_var_i)
_result_var = 0.0;
for (const auto i : index_range(_summed_vars))
{
_result_var.Add(_scale_factors[i], *_summed_vars[i]);
}
}

#endif
3 changes: 1 addition & 2 deletions test/tests/mfem/submeshes/cut_closed_coil.i
Original file line number Diff line number Diff line change
Expand Up @@ -141,8 +141,7 @@ coil_conductivity = 1.0
[update_total_e_field]
type = MFEMSumAux
variable = e_field
first_source_variable = induced_e_field
second_source_variable = external_e_field
source_variables = 'induced_e_field external_e_field'
execute_on = TIMESTEP_END
execution_order_group = 3
[]
Expand Down
3 changes: 1 addition & 2 deletions test/tests/mfem/submeshes/hphi_magnetostatic.i
Original file line number Diff line number Diff line change
Expand Up @@ -136,8 +136,7 @@ vacuum_permeability = 1.0
[update_total_h_field]
type = MFEMSumAux
variable = vacuum_h_field
first_source_variable = background_h_field
second_source_variable = cut_function_field
source_variables = 'background_h_field cut_function_field'
execute_on = TIMESTEP_END
execution_order_group = 3
[]
Expand Down
86 changes: 86 additions & 0 deletions unit/src/MFEMAuxKernelTest.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,86 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "MFEMObjectUnitTest.h"
#include "MFEMSumAux.h"

class MFEMAuxKernelTest : public MFEMObjectUnitTest
{
public:
MFEMAuxKernelTest() : MFEMObjectUnitTest("MooseUnitApp") {}
};

/**
* Test MFEMSumAux creates sums input GridFunctions successfully.
*/
TEST_F(MFEMAuxKernelTest, MFEMSumAux)
{
// Register dummy (Par)GridFunctions to test MFEMSumAux against
auto pm = _mfem_mesh_ptr->getMFEMParMeshPtr().get();
mfem::common::H1_ParFESpace fe(pm, 1);
mfem::common::H1_ParFESpace different_fe(pm, 2);

auto pgf_1 = std::make_shared<mfem::ParGridFunction>(&fe);
auto pgf_2 = std::make_shared<mfem::ParGridFunction>(&fe);
auto pgf_3 = std::make_shared<mfem::ParGridFunction>(&fe);
auto pgf_out = std::make_shared<mfem::ParGridFunction>(&fe);
auto pgf_ho = std::make_shared<mfem::ParGridFunction>(&different_fe);

_mfem_problem->getProblemData().gridfunctions.Register("source_variable_1", pgf_1);
_mfem_problem->getProblemData().gridfunctions.Register("source_variable_2", pgf_2);
_mfem_problem->getProblemData().gridfunctions.Register("source_variable_3", pgf_3);
_mfem_problem->getProblemData().gridfunctions.Register("summed_variable", pgf_out);
_mfem_problem->getProblemData().gridfunctions.Register("source_ho_variable", pgf_ho);

// Initialise variables that are to be summed over
mfem::ConstantCoefficient initial_condition_coef(2.0);
pgf_1->ProjectCoefficient(initial_condition_coef);
pgf_2->ProjectCoefficient(initial_condition_coef);
pgf_3->ProjectCoefficient(initial_condition_coef);
pgf_ho->ProjectCoefficient(initial_condition_coef);

{
// Construct auxkernel
InputParameters auxkernel_params = _factory.getValidParams("MFEMSumAux");
auxkernel_params.set<AuxVariableName>("variable") = "summed_variable";
auxkernel_params.set<std::vector<VariableName>>("source_variables") = {
"source_variable_1", "source_variable_2", "source_variable_3"};
auxkernel_params.set<std::vector<mfem::real_t>>("scale_factors") = {1.0, 2.0, 5.0};

MFEMSumAux & auxkernel = addObject<MFEMSumAux>("MFEMSumAux", "auxkernel1", auxkernel_params);
auxkernel.execute();

// Check the value of the output gridfunction is the scaled sum of the components
ASSERT_EQ(pgf_out->GetData()[5], 16.0);
}
{
// Check for failure if a variable to be summed has a different FESpace from the parent
InputParameters auxkernel_params = _factory.getValidParams("MFEMSumAux");
auxkernel_params.set<AuxVariableName>("variable") = "summed_variable";
auxkernel_params.set<std::vector<VariableName>>("source_variables") = {
"source_variable_1", "source_ho_variable", "source_variable_3"};
auxkernel_params.set<std::vector<mfem::real_t>>("scale_factors") = {1.0, 2.0, 5.0};
EXPECT_THROW(addObject<MFEMSumAux>("MFEMSumAux", "failed_auxkernel", auxkernel_params),
std::runtime_error);
}
{
// Check for failure if an inconsistent number of scale factors are provided
InputParameters auxkernel_params = _factory.getValidParams("MFEMSumAux");
auxkernel_params.set<AuxVariableName>("variable") = "summed_variable";
auxkernel_params.set<std::vector<VariableName>>("source_variables") = {
"source_variable_1", "source_variable_2", "source_variable_3"};
auxkernel_params.set<std::vector<mfem::real_t>>("scale_factors") = {1.0, 2.0};
EXPECT_THROW(addObject<MFEMSumAux>("MFEMSumAux", "failed_scaled_auxkernel", auxkernel_params),
std::runtime_error);
}
}

#endif
9 changes: 9 additions & 0 deletions unit/src/MFEMEssentialBCTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "libmesh/ignore_warnings.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMFESpaceTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "gtest/gtest.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMIntegratedBCTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "libmesh/ignore_warnings.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMKernelTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "MFEMObjectUnitTest.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMMaterialTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "MFEMObjectUnitTest.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMMeshTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "gtest/gtest.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMPostprocessorTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "MFEMObjectUnitTest.h"
Expand Down
9 changes: 9 additions & 0 deletions unit/src/MFEMSolverTest.C
Original file line number Diff line number Diff line change
@@ -1,3 +1,12 @@
//* This file is part of the MOOSE framework
//* https://mooseframework.inl.gov
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#ifdef MOOSE_MFEM_ENABLED

#include "MFEMObjectUnitTest.h"
Expand Down