diff --git a/include/userobjects/MFEMGeneralUserObject.h b/include/userobjects/MFEMGeneralUserObject.h index 907b6944..aea72da1 100644 --- a/include/userobjects/MFEMGeneralUserObject.h +++ b/include/userobjects/MFEMGeneralUserObject.h @@ -2,6 +2,7 @@ // MOOSE includes #include "GeneralUserObject.h" +#include "mfem.hpp" // Forwards declaration. class MFEMProblem; @@ -20,6 +21,11 @@ class MFEMGeneralUserObject : public GeneralUserObject // FIXME: Should this be marked as const if it is returning a non-const reference? MFEMProblem & getMFEMProblem() const { return _mfem_problem; } + /// Returns references to coefficients stored in the MFEMProblem PropertiesManager. + mfem::Coefficient & getScalarProperty(const std::string & name); + mfem::VectorCoefficient & getVectorProperty(const std::string & name); + mfem::MatrixCoefficient & getMatrixProperty(const std::string & name); + void execute() override {} void initialize() override {} diff --git a/src/bcs/MFEMConvectiveHeatFluxBC.C b/src/bcs/MFEMConvectiveHeatFluxBC.C index 09c4586c..89eb0e45 100644 --- a/src/bcs/MFEMConvectiveHeatFluxBC.C +++ b/src/bcs/MFEMConvectiveHeatFluxBC.C @@ -22,10 +22,9 @@ MFEMConvectiveHeatFluxBC::validParams() // TODO: Currently assumes the vector function coefficient is 3D MFEMConvectiveHeatFluxBC::MFEMConvectiveHeatFluxBC(const InputParameters & parameters) : MFEMIntegratedBC(parameters), - _heat_transfer_coef(getMFEMProblem().getProperties().getScalarProperty( + _heat_transfer_coef(getScalarProperty( getParam("heat_transfer_coefficient"))), - _T_inf_coef(getMFEMProblem().getProperties().getScalarProperty( - getParam("T_infinity"))), + _T_inf_coef(getScalarProperty(getParam("T_infinity"))), _external_heat_flux_coef(getMFEMProblem().makeScalarCoefficient( _heat_transfer_coef, _T_inf_coef)) { diff --git a/src/bcs/MFEMScalarFunctorDirichletBC.C b/src/bcs/MFEMScalarFunctorDirichletBC.C index 9a0a88c0..c9a0c671 100644 --- a/src/bcs/MFEMScalarFunctorDirichletBC.C +++ b/src/bcs/MFEMScalarFunctorDirichletBC.C @@ -14,7 +14,7 @@ MFEMScalarFunctorDirichletBC::validParams() MFEMScalarFunctorDirichletBC::MFEMScalarFunctorDirichletBC(const InputParameters & parameters) : MFEMEssentialBC(parameters), _coef_name(getParam("coefficient")), - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/bcs/MFEMVectorFunctorBoundaryIntegratedBC.C b/src/bcs/MFEMVectorFunctorBoundaryIntegratedBC.C index b44e9e5b..499f0849 100644 --- a/src/bcs/MFEMVectorFunctorBoundaryIntegratedBC.C +++ b/src/bcs/MFEMVectorFunctorBoundaryIntegratedBC.C @@ -15,8 +15,8 @@ MFEMVectorFunctorBoundaryIntegratedBC::validParams() MFEMVectorFunctorBoundaryIntegratedBC::MFEMVectorFunctorBoundaryIntegratedBC( const InputParameters & parameters) : MFEMIntegratedBC(parameters), - _vec_coef(getMFEMProblem().getProperties().getVectorProperty( - getParam("vector_coefficient"))) + _vec_coef( + getVectorProperty(getParam("vector_coefficient"))) { } diff --git a/src/bcs/MFEMVectorFunctorDirichletBCBase.C b/src/bcs/MFEMVectorFunctorDirichletBCBase.C index d4dcde86..4ab70188 100644 --- a/src/bcs/MFEMVectorFunctorDirichletBCBase.C +++ b/src/bcs/MFEMVectorFunctorDirichletBCBase.C @@ -14,7 +14,7 @@ MFEMVectorFunctorDirichletBCBase::validParams() MFEMVectorFunctorDirichletBCBase::MFEMVectorFunctorDirichletBCBase( const InputParameters & parameters) : MFEMEssentialBC(parameters), - _vec_coef(getMFEMProblem().getProperties().getVectorProperty( - getParam("vector_coefficient"))) + _vec_coef( + getVectorProperty(getParam("vector_coefficient"))) { } diff --git a/src/bcs/MFEMVectorFunctorNormalIntegratedBC.C b/src/bcs/MFEMVectorFunctorNormalIntegratedBC.C index b8b03e4c..8b9259e0 100644 --- a/src/bcs/MFEMVectorFunctorNormalIntegratedBC.C +++ b/src/bcs/MFEMVectorFunctorNormalIntegratedBC.C @@ -16,8 +16,8 @@ MFEMVectorFunctorNormalIntegratedBC::validParams() MFEMVectorFunctorNormalIntegratedBC::MFEMVectorFunctorNormalIntegratedBC( const InputParameters & parameters) : MFEMIntegratedBC(parameters), - _vec_coef(getMFEMProblem().getProperties().getVectorProperty( - getParam("vector_coefficient"))) + _vec_coef( + getVectorProperty(getParam("vector_coefficient"))) { } diff --git a/src/kernels/MFEMCurlCurlKernel.C b/src/kernels/MFEMCurlCurlKernel.C index 5aa4f00b..3b5853e1 100644 --- a/src/kernels/MFEMCurlCurlKernel.C +++ b/src/kernels/MFEMCurlCurlKernel.C @@ -22,7 +22,7 @@ MFEMCurlCurlKernel::MFEMCurlCurlKernel(const InputParameters & parameters) _coef_name(getParam("coefficient")), // FIXME: The MFEM bilinear form can also handle vector and matrix // coefficients, so ideally we'd handle all three too. - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMDiffusionKernel.C b/src/kernels/MFEMDiffusionKernel.C index 6b0e98d9..3a7c4b44 100644 --- a/src/kernels/MFEMDiffusionKernel.C +++ b/src/kernels/MFEMDiffusionKernel.C @@ -21,7 +21,7 @@ MFEMDiffusionKernel::MFEMDiffusionKernel(const InputParameters & parameters) _coef_name(getParam("coefficient")), // FIXME: The MFEM bilinear form can also handle vector and matrix // coefficients, so ideally we'd handle all three too. - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMDivDivKernel.C b/src/kernels/MFEMDivDivKernel.C index 37db579c..a1aedd01 100644 --- a/src/kernels/MFEMDivDivKernel.C +++ b/src/kernels/MFEMDivDivKernel.C @@ -24,7 +24,7 @@ MFEMDivDivKernel::MFEMDivDivKernel(const InputParameters & parameters) _coef_name(getParam("coefficient")), // FIXME: The MFEM bilinear form can also handle vector and matrix // coefficients, so ideally we'd handle all three too. - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMLinearElasticityKernel.C b/src/kernels/MFEMLinearElasticityKernel.C index 43873efc..72cf5b2c 100644 --- a/src/kernels/MFEMLinearElasticityKernel.C +++ b/src/kernels/MFEMLinearElasticityKernel.C @@ -28,8 +28,8 @@ MFEMLinearElasticityKernel::MFEMLinearElasticityKernel(const InputParameters & p : MFEMKernel(parameters), _lambda_name(getParam("lambda")), _mu_name(getParam("mu")), - _lambda(getMFEMProblem().getProperties().getScalarProperty(_lambda_name)), - _mu(getMFEMProblem().getProperties().getScalarProperty(_mu_name)) + _lambda(getScalarProperty(_lambda_name)), + _mu(getScalarProperty(_mu_name)) { } diff --git a/src/kernels/MFEMMassKernel.C b/src/kernels/MFEMMassKernel.C index 9ff2db01..26ac088b 100644 --- a/src/kernels/MFEMMassKernel.C +++ b/src/kernels/MFEMMassKernel.C @@ -18,7 +18,7 @@ MFEMMassKernel::validParams() MFEMMassKernel::MFEMMassKernel(const InputParameters & parameters) : MFEMKernel(parameters), _coef_name(getParam("coefficient")), - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMMixedVectorGradientKernel.C b/src/kernels/MFEMMixedVectorGradientKernel.C index f340216a..870f905b 100644 --- a/src/kernels/MFEMMixedVectorGradientKernel.C +++ b/src/kernels/MFEMMixedVectorGradientKernel.C @@ -20,7 +20,7 @@ MFEMMixedVectorGradientKernel::MFEMMixedVectorGradientKernel(const InputParamete _coef_name(getParam("coefficient")), // FIXME: The MFEM bilinear form can also handle vector and matrix // coefficients, so ideally we'd handle all three too. - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMVectorFEDomainLFKernel.C b/src/kernels/MFEMVectorFEDomainLFKernel.C index 9f54be75..46469d23 100644 --- a/src/kernels/MFEMVectorFEDomainLFKernel.C +++ b/src/kernels/MFEMVectorFEDomainLFKernel.C @@ -17,8 +17,8 @@ MFEMVectorFEDomainLFKernel::validParams() MFEMVectorFEDomainLFKernel::MFEMVectorFEDomainLFKernel(const InputParameters & parameters) : MFEMKernel(parameters), - _vec_coef(getMFEMProblem().getProperties().getVectorProperty( - getParam("vector_coefficient"))) + _vec_coef( + getVectorProperty(getParam("vector_coefficient"))) { } diff --git a/src/kernels/MFEMVectorFEMassKernel.C b/src/kernels/MFEMVectorFEMassKernel.C index 74cbced7..611f5893 100644 --- a/src/kernels/MFEMVectorFEMassKernel.C +++ b/src/kernels/MFEMVectorFEMassKernel.C @@ -20,7 +20,7 @@ MFEMVectorFEMassKernel::MFEMVectorFEMassKernel(const InputParameters & parameter _coef_name(getParam("coefficient")), // FIXME: The MFEM bilinear form can also handle vector and matrix // coefficients, so ideally we'd handle all three too. - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/kernels/MFEMVectorFEWeakDivergenceKernel.C b/src/kernels/MFEMVectorFEWeakDivergenceKernel.C index 3a67b45e..ca7fea7a 100644 --- a/src/kernels/MFEMVectorFEWeakDivergenceKernel.C +++ b/src/kernels/MFEMVectorFEWeakDivergenceKernel.C @@ -19,7 +19,7 @@ MFEMVectorFEWeakDivergenceKernel::MFEMVectorFEWeakDivergenceKernel( const InputParameters & parameters) : MFEMKernel(parameters), _coef_name(getParam("coefficient")), - _coef(getMFEMProblem().getProperties().getScalarProperty(_coef_name)) + _coef(getScalarProperty(_coef_name)) { } diff --git a/src/userobjects/MFEMGeneralUserObject.C b/src/userobjects/MFEMGeneralUserObject.C index 5ba78179..3bdab4b4 100644 --- a/src/userobjects/MFEMGeneralUserObject.C +++ b/src/userobjects/MFEMGeneralUserObject.C @@ -15,3 +15,21 @@ MFEMGeneralUserObject::MFEMGeneralUserObject(const InputParameters & parameters) : GeneralUserObject(parameters), _mfem_problem(static_cast(_fe_problem)) { } + +mfem::Coefficient & +MFEMGeneralUserObject::getScalarProperty(const std::string & name) +{ + return getMFEMProblem().getProperties().getScalarProperty(name); +} + +mfem::VectorCoefficient & +MFEMGeneralUserObject::getVectorProperty(const std::string & name) +{ + return getMFEMProblem().getProperties().getVectorProperty(name); +} + +mfem::MatrixCoefficient & +MFEMGeneralUserObject::getMatrixProperty(const std::string & name) +{ + return getMFEMProblem().getProperties().getMatrixProperty(name); +}