Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\/ M anipulation |
-------------------------------------------------------------------------------

Kutalmış Berçin
committed
Copyright (C) 2019-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "nutUWallFunctionFvPatchScalarField.H"
#include "turbulenceModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::scalarField>
Foam::nutUWallFunctionFvPatchScalarField::calcNut() const
{
const label patchi = patch().index();
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);

Andrew Heather
committed
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
const tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
tmp<scalarField> tyPlus = calcYPlus(magUp);
tmp<scalarField> tnutw(new scalarField(patch().size(), Zero));
scalarField& nutw = tnutw.ref();
forAll(yPlus, facei)
{
// Viscous sublayer contribution
const scalar nutVis = 0.0;
// Inertial sublayer contribution
const scalar nutLog =

Kutalmış Berçin
committed
nuw[facei]
*(yPlus[facei]*kappa_/log(max(E_*yPlus[facei], 1 + 1e-4)) - 1.0);
nutw[facei] = blend(nutVis, nutLog, yPlus[facei]);
}
return tnutw;
}
Foam::tmp<Foam::scalarField>
Foam::nutUWallFunctionFvPatchScalarField::calcYPlus
(
const scalarField& magUp
) const
{
const label patchi = patch().index();
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
const scalarField& y = turbModel.y()[patchi];
const tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
tmp<scalarField> tyPlus(new scalarField(patch().size(), Zero));
scalarField& yPlus = tyPlus.ref();
forAll(yPlus, facei)
{
const scalar kappaRe = kappa_*magUp[facei]*y[facei]/nuw[facei];
scalar yp = yPlusLam_;
const scalar ryPlusLam = 1.0/yp;
int iter = 0;
scalar yPlusLast = 0.0;
do
{
yPlusLast = yp;
yp = (kappaRe + yp)/(1.0 + log(E_*yp));
} while (mag(ryPlusLam*(yp - yPlusLast)) > 0.01 && ++iter < 10 );
yPlus[facei] = max(0.0, yp);
}
return tyPlus;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
nutWallFunctionFvPatchScalarField(p, iF)
{}
Foam::nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField
(
const nutUWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
{}
Foam::nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
nutWallFunctionFvPatchScalarField(p, iF, dict)
{}
Foam::nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField
(
const nutUWallFunctionFvPatchScalarField& sawfpsf
)
:
nutWallFunctionFvPatchScalarField(sawfpsf)
{}
Foam::nutUWallFunctionFvPatchScalarField::nutUWallFunctionFvPatchScalarField
(
const nutUWallFunctionFvPatchScalarField& sawfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
nutWallFunctionFvPatchScalarField(sawfpsf, iF)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::scalarField>
Foam::nutUWallFunctionFvPatchScalarField::yPlus() const
{
const label patchi = patch().index();

Kutalmış Berçin
committed
const auto& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);

Kutalmış Berçin
committed
const scalarField& y = turbModel.y()[patchi];
tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
tmp<scalarField> tnuEff = turbModel.nuEff(patchi);
const scalarField& nuEff = tnuEff();

Andrew Heather
committed
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
const scalarField magUp(mag(Uw.patchInternalField() - Uw));

Kutalmış Berçin
committed
const scalarField magGradUw(mag(Uw.snGrad()));

Kutalmış Berçin
committed
tmp<scalarField> tyPlus = calcYPlus(magUp);
scalarField& yPlus = tyPlus.ref();
forAll(yPlus, facei)
{
if (yPlusLam_ > yPlus[facei])
{
// viscous sublayer
yPlus[facei] =
y[facei]*sqrt(nuEff[facei]*magGradUw[facei])/nuw[facei];
}
}
return tyPlus;
}
void Foam::nutUWallFunctionFvPatchScalarField::write
(
Ostream& os
) const
{
fvPatchField<scalar>::write(os);
writeLocalEntries(os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
nutUWallFunctionFvPatchScalarField
);
}
// ************************************************************************* //