ENH: Using specie::Tstd instead of explicitly declaring (potentially different) Tstd

This commit is contained in:
andy
2010-09-07 16:10:16 +01:00
parent 5e9c1db849
commit 0de24cd710
2 changed files with 12 additions and 6 deletions

View File

@ -30,6 +30,7 @@ License
#include "addToRunTimeSelectionTable.H" #include "addToRunTimeSelectionTable.H"
#include "zeroGradientFvPatchFields.H" #include "zeroGradientFvPatchFields.H"
#include "directMappedWallPolyPatch.H" #include "directMappedWallPolyPatch.H"
#include "specie.H"
// Sub-models // Sub-models
#include "heatTransferModel.H" #include "heatTransferModel.H"
@ -195,12 +196,12 @@ Foam::tmp<Foam::fvScalarMatrix> Foam::surfaceFilmModels::thermoSingleLayer::q
volScalarField& hs volScalarField& hs
) const ) const
{ {
const dimensionedScalar Tstd("Tstd", dimTemperature, specie::Tstd);
return return
( (
- fvm::Sp(htcs_->h()/cp_, hs) - fvm::Sp(htcs_->h()/cp_, hs) - htcs_->h()*(Tstd - Ts_)
- htcs_->h()*(dimensionedScalar("Tstd", dimTemperature, 298.15) - Ts_) - fvm::Sp(htcw_->h()/cp_, hs) - htcw_->h()*(Tstd - Tw_)
- fvm::Sp(htcw_->h()/cp_, hs)
- htcw_->h()*(dimensionedScalar("Tstd", dimTemperature, 298.15) - Tw_)
); );
} }

View File

@ -24,6 +24,7 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "thermoSingleLayer.H" #include "thermoSingleLayer.H"
#include "specie.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -46,6 +47,8 @@ Foam::surfaceFilmModels::thermoSingleLayer::hs
const volScalarField& T const volScalarField& T
) const ) const
{ {
const dimensionedScalar Tstd("Tstd", dimTemperature, specie::Tstd);
return tmp<volScalarField> return tmp<volScalarField>
( (
new volScalarField new volScalarField
@ -58,7 +61,7 @@ Foam::surfaceFilmModels::thermoSingleLayer::hs
IOobject::NO_READ, IOobject::NO_READ,
IOobject::NO_WRITE IOobject::NO_WRITE
), ),
cp_*(T - (dimensionedScalar("Tstd", dimTemperature, 298.15))), cp_*(T - Tstd),
zeroGradientFvPatchScalarField::typeName zeroGradientFvPatchScalarField::typeName
) )
); );
@ -71,6 +74,8 @@ Foam::surfaceFilmModels::thermoSingleLayer::T
const volScalarField& hs const volScalarField& hs
) const ) const
{ {
const dimensionedScalar Tstd("Tstd", dimTemperature, specie::Tstd);
return tmp<volScalarField> return tmp<volScalarField>
( (
new volScalarField new volScalarField
@ -83,7 +88,7 @@ Foam::surfaceFilmModels::thermoSingleLayer::T
IOobject::NO_READ, IOobject::NO_READ,
IOobject::NO_WRITE IOobject::NO_WRITE
), ),
hs/cp_ + dimensionedScalar("Tstd", dimTemperature, 298.15), hs/cp_ + Tstd,
zeroGradientFvPatchScalarField::typeName zeroGradientFvPatchScalarField::typeName
) )
); );