Files
OpenFOAM-12/applications/modules/multiphaseEuler/phaseSystems/OneResistanceHeatTransferPhaseSystem/OneResistanceHeatTransferPhaseSystem.C
Will Bainbridge 597121a4a7 multiphaseEuler: Library reorganisation
This change makes multiphaseEuler more consistent with other modules and
makes its sub-libraries less inter-dependent. Some left-over references
to multiphaseEulerFoam have also been removed.
2023-09-15 14:45:26 +01:00

134 lines
3.8 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2020-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 "OneResistanceHeatTransferPhaseSystem.H"
#include "heatTransferModel.H"
#include "fvmSup.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::OneResistanceHeatTransferPhaseSystem<BasePhaseSystem>::
OneResistanceHeatTransferPhaseSystem
(
const fvMesh& mesh
)
:
HeatTransferPhaseSystem<BasePhaseSystem>(mesh)
{
this->generateInterfacialModels(heatTransferModels_);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::OneResistanceHeatTransferPhaseSystem<BasePhaseSystem>::
~OneResistanceHeatTransferPhaseSystem()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::autoPtr<Foam::phaseSystem::heatTransferTable>
Foam::OneResistanceHeatTransferPhaseSystem<BasePhaseSystem>::
heatTransfer() const
{
autoPtr<phaseSystem::heatTransferTable> eqnsPtr
(
new phaseSystem::heatTransferTable()
);
phaseSystem::heatTransferTable& eqns = eqnsPtr();
forAll(this->phaseModels_, phasei)
{
const phaseModel& phase = this->phaseModels_[phasei];
eqns.insert
(
phase.name(),
new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
);
}
// Heat transfer across the interface
forAllConstIter
(
heatTransferModelTable,
heatTransferModels_,
heatTransferModelIter
)
{
const phaseInterface interface(*this, heatTransferModelIter.key());
const volScalarField K(heatTransferModelIter()->K());
forAllConstIter(phaseInterface, interface, iter)
{
const phaseModel& phase = iter();
const phaseModel& otherPhase = iter.otherPhase();
const volScalarField& he = phase.thermo().he();
const volScalarField Cpv(phase.thermo().Cpv());
const volScalarField phaseK
(
iter.otherPhase()
/max(iter.otherPhase(), iter.otherPhase().residualAlpha())
*K
);
*eqns[phase.name()] +=
phaseK*(otherPhase.thermo().T() - phase.thermo().T() + he/Cpv)
- fvm::Sp(phaseK/Cpv, he);
}
}
return eqnsPtr;
}
template<class BasePhaseSystem>
bool Foam::OneResistanceHeatTransferPhaseSystem<BasePhaseSystem>::read()
{
if (HeatTransferPhaseSystem<BasePhaseSystem>::read())
{
bool readOK = true;
// Models ...
return readOK;
}
else
{
return false;
}
}
// ************************************************************************* //