Skip to content
Snippets Groups Projects
Commit 565825b7 authored by andy's avatar andy
Browse files

ENH: Surface film - heat transfer only active where film is present

parent 010e2e2d
Branches
Tags
No related merge requests found
......@@ -252,10 +252,22 @@ tmp<fvScalarMatrix> thermoSingleLayer::q(volScalarField& hs) const
{
dimensionedScalar Tstd("Tstd", dimTemperature, 298.15);
volScalarField htcst(htcs_->h());
volScalarField htcwt(htcw_->h());
forAll(alpha_, i)
{
htcst[i] *= max(alpha_[i], ROOTVSMALL);
htcwt[i] *= max(alpha_[i], ROOTVSMALL);
}
htcst.correctBoundaryConditions();
htcwt.correctBoundaryConditions();
return
(
- fvm::Sp(htcs_->h()/Cp_, hs) - htcs_->h()*(Tstd - TPrimary_)
- fvm::Sp(htcw_->h()/Cp_, hs) - htcw_->h()*(Tstd - Tw_)
- fvm::Sp(htcst/Cp_, hs) - htcst*(Tstd - TPrimary_)
- fvm::Sp(htcwt/Cp_, hs) - htcwt*(Tstd - Tw_)
);
}
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment