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;
        }
    }
This commit is contained in:
Will Bainbridge
2023-11-28 12:33:21 +00:00
parent 3a6f9029d3
commit 01d0af39be
3 changed files with 214 additions and 0 deletions

View File

@ -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

View File

@ -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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#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::volScalarField>
Foam::heatTransferModels::Prandtl::K
(
const scalar residualAlpha
) const
{
const dragModel& drag =
interface_.fluid().lookupInterfacialModel<dragModel>(interface_);
const volScalarField CpAvg
(
interface_.thermo1().Cp()*interface_.thermo2().Cp()
/(interface_.thermo1().Cp() + interface_.thermo2().Cp())
);
return drag.K()*CpAvg/Pr_;
}
// ************************************************************************* //

View File

@ -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 <http://www.gnu.org/licenses/>.
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<phaseInterface> 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<volScalarField> K(const scalar residualAlpha) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace heatTransferModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //