Skip to content
Snippets Groups Projects
Commit 06e4f4db authored by Henry's avatar Henry
Browse files

thermoSingleLayer: revert change to q function made shortly before the release of OpenFOAM-2.3.1

Resolves bug-report http://www.openfoam.org/mantisbt/view.php?id=1502
parent 52eb0c83
Branches
Tags
No related merge requests found
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
......@@ -257,16 +257,9 @@ tmp<fvScalarMatrix> thermoSingleLayer::q(volScalarField& hs) const
{
dimensionedScalar Tstd("Tstd", dimTemperature, 298.15);
volScalarField htcst(htcs_->h());
volScalarField htcwt(htcw_->h());
const volScalarField mask(pos(delta_ - deltaSmall_));
forAll(mask, i)
{
htcst[i] *= max(mask[i], ROOTVSMALL);
htcwt[i] *= max(mask[i], ROOTVSMALL);
}
volScalarField boundedAlpha(max(alpha_, ROOTVSMALL));
volScalarField htcst(htcs_->h()*boundedAlpha);
volScalarField htcwt(htcw_->h()*boundedAlpha);
htcst.correctBoundaryConditions();
htcwt.correctBoundaryConditions();
......
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