Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2016-2023 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 "yPlus.H"
#include "volFields.H"
#include "turbulenceModel.H"
#include "nutWallFunctionFvPatchScalarField.H"
#include "wallFvPatch.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
Henry Weller
committed
{
namespace functionObjects
defineTypeNameAndDebug(yPlus, 0);
addToRunTimeSelectionTable(functionObject, yPlus, dictionary);
Henry Weller
committed
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::functionObjects::yPlus::writeFileHeader(Ostream& os) const
writeHeader(os, "y+ ()");
writeCommented(os, "Time");
writeTabbed(os, "patch");
writeTabbed(os, "min");
writeTabbed(os, "max");
writeTabbed(os, "average");
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Henry Weller
committed
Foam::functionObjects::yPlus::yPlus
const Time& runTime,
const dictionary& dict
fvMeshFunctionObject(name, runTime, dict),
writeFile(obr_, name, typeName, dict),
useWallFunction_(true),
writeFields_(true) // May change in the future
writeFileHeader(file());
Henry Weller
committed
(
new volScalarField
Henry Weller
committed
IOobject
scopedName(typeName),
mesh_.time().timeName(),
mesh_,
Henry Weller
committed
IOobject::NO_READ,
IOobject::NO_WRITE,
IOobject::REGISTER
Henry Weller
committed
),
mesh_,
dimensionedScalar(dimless, Zero)
Henry Weller
committed
)
);
Henry Weller
committed
}
Henry Weller
committed
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::functionObjects::yPlus::read(const dictionary& dict)
if (fvMeshFunctionObject::read(dict) && writeFile::read(dict))
{
useWallFunction_ = true;
writeFields_ = true; // May change in the future
dict.readIfPresent("useWallFunction", useWallFunction_);
dict.readIfPresent("writeFields", writeFields_);
return true;
}
return false;
bool Foam::functionObjects::yPlus::execute()
auto& yPlus = lookupObjectRef<volScalarField>(scopedName(typeName));
if (foundObject<turbulenceModel>(turbulenceModel::propertiesName))
volScalarField::Boundary& yPlusBf = yPlus.boundaryFieldRef();
const turbulenceModel& model =
lookupObject<turbulenceModel>
(
turbulenceModel::propertiesName
);
Henry Weller
committed
const nearWallDist nwd(mesh_);
const volScalarField::Boundary& d = nwd.y();
// nut needed for wall function patches
tmp<volScalarField> tnut = model.nut();
const volScalarField::Boundary& nutBf = tnut().boundaryField();
const volVectorField::Boundary& UBf = model.U().boundaryField();
const fvPatchList& patches = mesh_.boundary();
forAll(patches, patchi)
{
const fvPatch& patch = patches[patchi];
const auto* nutWallPatch =
isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]);
if (useWallFunction_ && nutWallPatch)
{
yPlusBf[patchi] = nutWallPatch->yPlus();
}
else if (isA<wallFvPatch>(patch))
{
yPlusBf[patchi] =
d[patchi]
*sqrt
(
*mag(UBf[patchi].snGrad())
}
}
Henry Weller
committed
}
else
{
Henry Weller
committed
<< "Unable to find turbulence model in the "
<< "database: yPlus will not be calculated" << endl;

Kutalmış Berçin
committed
if (postProcess)
{
WarningInFunction
<< "Please try to use the solver option -postProcess, e.g.:"
<< " <solver> -postProcess -func yPlus" << endl;
}
return true;
bool Foam::functionObjects::yPlus::write()
const auto& yPlus = obr_.lookupObject<volScalarField>(scopedName(typeName));
Log << type() << ' ' << name() << " write:" << nl;
if (writeFields_)
{
Log << " writing field " << yPlus.name() << endl;
yPlus.write();
}
Henry Weller
committed
const volScalarField::Boundary& yPlusBf = yPlus.boundaryField();
const fvPatchList& patches = mesh_.boundary();
Henry Weller
committed
forAll(patches, patchi)
{
const fvPatch& patch = patches[patchi];
if (isA<wallFvPatch>(patch))
{
const scalarField& yPlusp = yPlusBf[patchi];
const scalar minYplus = gMin(yPlusp);
const scalar maxYplus = gMax(yPlusp);
const scalar avgYplus = gAverage(yPlusp);
if (UPstream::master())
Henry Weller
committed
{

Andrew Heather
committed
writeCurrentTime(file());
Henry Weller
committed
file()
<< token::TAB << patch.name()
Henry Weller
committed
<< token::TAB << minYplus
<< token::TAB << maxYplus
<< token::TAB << avgYplus
<< endl;
}
Log << " patch " << patch.name()
<< " y+ : min = " << minYplus
<< ", max = " << maxYplus
<< ", average = " << avgYplus << endl;
Henry Weller
committed
}
}
return true;