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 "MFEMConvectiveHeatFluxBC.h" 13 : #include "MFEMProblem.h" 14 : 15 : registerMooseObject("MooseApp", MFEMConvectiveHeatFluxBC); 16 : 17 : InputParameters 18 8648 : MFEMConvectiveHeatFluxBC::validParams() 19 : { 20 8648 : InputParameters params = MFEMIntegratedBC::validParams(); 21 : // FIXME: Should these really be specified via properties? T_infinity in particular? Use functions 22 : // instead? 23 8648 : params.addClassDescription( 24 : "Convective heat transfer boundary condition with temperature and heat " 25 : "transfer coefficent given by material properties to add to MFEM problems."); 26 8648 : params.addParam<MFEMScalarCoefficientName>( 27 : "T_infinity", "0.", "Name of a coefficient specifying the far-field temperature"); 28 8648 : params.addParam<MFEMScalarCoefficientName>( 29 : "heat_transfer_coefficient", 30 : "1.", 31 : "Name of the coefficient specifying the heat transfer coefficient"); 32 8648 : return params; 33 0 : } 34 : 35 9 : MFEMConvectiveHeatFluxBC::MFEMConvectiveHeatFluxBC(const InputParameters & parameters) 36 : : MFEMIntegratedBC(parameters), 37 9 : _heat_transfer_coef(getScalarCoefficient("heat_transfer_coefficient")), 38 9 : _T_inf_coef(getScalarCoefficient("T_infinity")), 39 9 : _external_heat_flux_coef( 40 18 : getMFEMProblem().getCoefficients().declareScalar<mfem::ProductCoefficient>( 41 18 : "__ConvectiveHeatFluxBC_" + parameters.get<std::string>("_unique_name"), 42 : _heat_transfer_coef, 43 9 : _T_inf_coef)) 44 : { 45 9 : } 46 : 47 : // Create a new MFEM integrator to apply to the RHS of the weak form. Ownership managed by the 48 : // caller. 49 : mfem::LinearFormIntegrator * 50 33 : MFEMConvectiveHeatFluxBC::createLFIntegrator() 51 : { 52 33 : return new mfem::BoundaryLFIntegrator(_external_heat_flux_coef); 53 : } 54 : 55 : // Create a new MFEM integrator to apply to LHS of the weak form. Ownership managed by the caller. 56 : mfem::BilinearFormIntegrator * 57 65 : MFEMConvectiveHeatFluxBC::createBFIntegrator() 58 : { 59 65 : return new mfem::BoundaryMassIntegrator(_heat_transfer_coef); 60 : } 61 : 62 : #endif