mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Added on-the-fly heat flux calculations to thermoFilm
This commit is contained in:
@ -370,6 +370,15 @@ public:
|
||||
inline const filmRadiationModel& radiation() const;
|
||||
|
||||
|
||||
// Derived fields (calculated on-the-fly)
|
||||
|
||||
//- Return the convective heat energy from film to wall
|
||||
inline tmp<scalarField> Qconvw(const label patchI) const;
|
||||
|
||||
//- Return the convective heat energy from primary region to film
|
||||
inline tmp<scalarField> Qconvp(const label patchI) const;
|
||||
|
||||
|
||||
// Evolution
|
||||
|
||||
//- Pre-evolve film hook
|
||||
|
||||
@ -24,6 +24,7 @@ License
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "thermoSingleLayer.H"
|
||||
#include "heatTransferModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -160,6 +161,26 @@ inline const filmRadiationModel& thermoSingleLayer::radiation() const
|
||||
}
|
||||
|
||||
|
||||
inline tmp<scalarField> thermoSingleLayer::Qconvw(const label patchI) const
|
||||
{
|
||||
const scalarField htc(htcw_->h()().boundaryField()[patchI]);
|
||||
const scalarField& Tp = T_.boundaryField()[patchI];
|
||||
const scalarField& Twp = Tw_.boundaryField()[patchI];
|
||||
|
||||
return htc*(Tp - Twp);
|
||||
}
|
||||
|
||||
|
||||
inline tmp<scalarField> thermoSingleLayer::Qconvp(const label patchI) const
|
||||
{
|
||||
const scalarField htc(htcs_->h()().boundaryField()[patchI]);
|
||||
const scalarField& Tp = T_.boundaryField()[patchI];
|
||||
const scalarField& Tpp = TPrimary_.boundaryField()[patchI];
|
||||
|
||||
return htc*(Tp - Tpp);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace surfaceFilmModels
|
||||
|
||||
Reference in New Issue
Block a user