From 01d0af39be1ea4b4b1b440dd7d6f220ecbef3ff9 Mon Sep 17 00:00:00 2001 From: Will Bainbridge Date: Tue, 28 Nov 2023 12:33:21 +0000 Subject: [PATCH] multiphaseEuler: Added Prandtl heat transfer model This simple model creates a heat transfer coefficient in proportion with the corresponding drag model's momentum transfer coefficient. A user-defined Prandtl number and a harmonic average of the phases' specific heats are used to specify the constant of proportionality. This model has no physical basis. It exists primarily for testing purposes. It has the advantage of being applicable to any interface, including those representing segregated configurations. Example usage: heatTransfer { gas_segregatedWith_liquid { type Prandtl; Pr 0.7; } } --- .../interfacialModels/Make/files | 1 + .../heatTransferModels/Prandtl/Prandtl.C | 90 +++++++++++++ .../heatTransferModels/Prandtl/Prandtl.H | 123 ++++++++++++++++++ 3 files changed, 214 insertions(+) create mode 100644 applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.C create mode 100644 applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.H diff --git a/applications/modules/multiphaseEuler/interfacialModels/Make/files b/applications/modules/multiphaseEuler/interfacialModels/Make/files index 0b3a539f95..d7774cd201 100644 --- a/applications/modules/multiphaseEuler/interfacialModels/Make/files +++ b/applications/modules/multiphaseEuler/interfacialModels/Make/files @@ -44,6 +44,7 @@ heatTransferModels/Gunn/Gunn.C heatTransferModels/sphericalHeatTransfer/sphericalHeatTransfer.C heatTransferModels/nonSphericalHeatTransfer/nonSphericalHeatTransfer.C heatTransferModels/timeScaleFilteredHeatTransfer/timeScaleFilteredHeatTransfer.C +heatTransferModels/Prandtl/Prandtl.C virtualMassModels/virtualMassModel/virtualMassModel.C virtualMassModels/virtualMassModel/virtualMassModelNew.C diff --git a/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.C b/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.C new file mode 100644 index 0000000000..052ebbbbbd --- /dev/null +++ b/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.C @@ -0,0 +1,90 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | Website: https://openfoam.org + \\ / A nd | Copyright (C) 2023 OpenFOAM Foundation + \\/ M anipulation | +------------------------------------------------------------------------------- +License + This file is part of OpenFOAM. + + OpenFOAM is free software: you can redistribute it and/or modify it + under the terms of the GNU General Public License as published by + the Free Software Foundation, either version 3 of the License, or + (at your option) any later version. + + OpenFOAM is distributed in the hope that it will be useful, but WITHOUT + ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or + FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License + for more details. + + You should have received a copy of the GNU General Public License + along with OpenFOAM. If not, see . + +\*---------------------------------------------------------------------------*/ + +#include "Prandtl.H" +#include "dragModel.H" +#include "addToRunTimeSelectionTable.H" + +// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // + +namespace Foam +{ +namespace heatTransferModels +{ + defineTypeNameAndDebug(Prandtl, 0); + addToRunTimeSelectionTable + ( + heatTransferModel, + Prandtl, + dictionary + ); +} +} + + +// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // + +Foam::heatTransferModels::Prandtl::Prandtl +( + const dictionary& dict, + const phaseInterface& interface, + const bool registerObject +) +: + heatTransferModel(dict, interface, registerObject), + interfacePtr_(interface.clone()), + interface_(interfacePtr_()), + Pr_("Pr", dimless, dict) +{} + + +// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // + +Foam::heatTransferModels::Prandtl::~Prandtl() +{} + + +// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // + +Foam::tmp +Foam::heatTransferModels::Prandtl::K +( + const scalar residualAlpha +) const +{ + const dragModel& drag = + interface_.fluid().lookupInterfacialModel(interface_); + + const volScalarField CpAvg + ( + interface_.thermo1().Cp()*interface_.thermo2().Cp() + /(interface_.thermo1().Cp() + interface_.thermo2().Cp()) + ); + + return drag.K()*CpAvg/Pr_; +} + + +// ************************************************************************* // diff --git a/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.H b/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.H new file mode 100644 index 0000000000..168f1339ae --- /dev/null +++ b/applications/modules/multiphaseEuler/interfacialModels/heatTransferModels/Prandtl/Prandtl.H @@ -0,0 +1,123 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | Website: https://openfoam.org + \\ / A nd | Copyright (C) 2023 OpenFOAM Foundation + \\/ M anipulation | +------------------------------------------------------------------------------- +License + This file is part of OpenFOAM. + + OpenFOAM is free software: you can redistribute it and/or modify it + under the terms of the GNU General Public License as published by + the Free Software Foundation, either version 3 of the License, or + (at your option) any later version. + + OpenFOAM is distributed in the hope that it will be useful, but WITHOUT + ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or + FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License + for more details. + + You should have received a copy of the GNU General Public License + along with OpenFOAM. If not, see . + +Class + Foam::heatTransferModels::Prandtl + +Description + This simple model creates a heat transfer coefficient in proportion with + the corresponding drag model's momentum transfer coefficient. A + user-defined Prandtl number and a harmonic average of the phases' + specific heats are used to specify the constant of proportionality. + + This model has no physical basis. It exists primarily for testing + purposes. It has the advantage of being applicable to any interface, + including those representing segregated configurations. + + Example usage: + \verbatim + heatTransfer + { + gas_segregatedWith_liquid + { + type Prandtl; + Pr 0.7; + } + } + \endverbatim + +SourceFiles + Prandtl.C + +\*---------------------------------------------------------------------------*/ + +#ifndef Prandtl_H +#define Prandtl_H + +#include "heatTransferModel.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ +namespace heatTransferModels +{ + +/*---------------------------------------------------------------------------*\ + Class Prandtl Declaration +\*---------------------------------------------------------------------------*/ + +class Prandtl +: + public heatTransferModel +{ + // Private Data + + //- Interface pointer + const autoPtr interfacePtr_; + + //- Interface + const phaseInterface& interface_; + + //- Prandtl number + const dimensionedScalar Pr_; + + +public: + + //- Runtime type information + TypeName("Prandtl"); + + + // Constructors + + //- Construct from a dictionary and an interface + Prandtl + ( + const dictionary& dict, + const phaseInterface& interface, + const bool registerObject + ); + + + //- Destructor + virtual ~Prandtl(); + + + // Member Functions + + //- The heat transfer function K used in the enthalpy equation + tmp K(const scalar residualAlpha) const; +}; + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace heatTransferModels +} // End namespace Foam + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +#endif + +// ************************************************************************* //