openfoam/applications/utilities/preProcessing/applyBoundaryLayer/applyBoundaryLayer.C
Henry 2aec249647 Updated the whole of OpenFOAM to use the new templated TurbulenceModels library
The old separate incompressible and compressible libraries have been removed.

Most of the commonly used RANS and LES models have been upgraded to the
new framework but there are a few missing which will be added over the
next few days, in particular the realizable k-epsilon model.  Some of
the less common incompressible RANS models have been introduced into the
new library instantiated for incompressible flow only.  If they prove to
be generally useful they can be templated for compressible and
multiphase application.

The Spalart-Allmaras DDES and IDDES models have been thoroughly
debugged, removing serious errors concerning the use of S rather than
Omega.

The compressible instances of the models have been augmented by a simple
backward-compatible eddyDiffusivity model for thermal transport based on
alphat and alphaEff.  This will be replaced with a separate run-time
selectable thermal transport model framework in a few weeks.

For simplicity and ease of maintenance and further development the
turbulent transport and wall modeling is based on nut/nuEff rather than
mut/muEff for compressible models so that all forms of turbulence models
can use the same wall-functions and other BCs.

All turbulence model selection made in the constant/turbulenceProperties
dictionary with RAS and LES as sub-dictionaries rather than in separate
files which added huge complexity for multiphase.

All tutorials have been updated so study the changes and update your own
cases by comparison with similar cases provided.

Sorry for the inconvenience in the break in backward-compatibility but
this update to the turbulence modeling is an essential step in the
future of OpenFOAM to allow more models to be added and maintained for a
wider range of cases and physics.  Over the next weeks and months more
turbulence models will be added of single and multiphase flow, more
additional sub-models and further development and testing of existing
models.  I hope this brings benefits to all OpenFOAM users.

Henry G. Weller
2015-01-21 19:21:39 +00:00

239 lines
7.0 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
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/>.
Application
applyBoundaryLayer
Description
Apply a simplified boundary-layer model to the velocity and
turbulence fields based on the 1/7th power-law.
The uniform boundary-layer thickness is either provided via the -ybl option
or calculated as the average of the distance to the wall scaled with
the thickness coefficient supplied via the option -Cbl. If both options
are provided -ybl is used.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "turbulentTransportModel.H"
#include "wallDist.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// turbulence constants - file-scope
static const scalar Cmu(0.09);
static const scalar kappa(0.41);
int main(int argc, char *argv[])
{
argList::addNote
(
"apply a simplified boundary-layer model to the velocity and\n"
"turbulence fields based on the 1/7th power-law."
);
argList::addOption
(
"ybl",
"scalar",
"specify the boundary-layer thickness"
);
argList::addOption
(
"Cbl",
"scalar",
"boundary-layer thickness as Cbl * mean distance to wall"
);
argList::addBoolOption
(
"writenut",
"write nut field"
);
#include "setRootCase.H"
if (!args.optionFound("ybl") && !args.optionFound("Cbl"))
{
FatalErrorIn(args.executable())
<< "Neither option 'ybl' or 'Cbl' have been provided to calculate "
<< "the boundary-layer thickness.\n"
<< "Please choose either 'ybl' OR 'Cbl'."
<< exit(FatalError);
}
else if (args.optionFound("ybl") && args.optionFound("Cbl"))
{
FatalErrorIn(args.executable())
<< "Both 'ybl' and 'Cbl' have been provided to calculate "
<< "the boundary-layer thickness.\n"
<< "Please choose either 'ybl' OR 'Cbl'."
<< exit(FatalError);
}
#include "createTime.H"
#include "createMesh.H"
#include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Modify velocity by applying a 1/7th power law boundary-layer
// u/U0 = (y/ybl)^(1/7)
// assumes U0 is the same as the current cell velocity
Info<< "Setting boundary layer velocity" << nl << endl;
scalar yblv = ybl.value();
forAll(U, cellI)
{
if (y[cellI] <= yblv)
{
mask[cellI] = 1;
U[cellI] *= ::pow(y[cellI]/yblv, (1.0/7.0));
}
}
mask.correctBoundaryConditions();
Info<< "Writing U\n" << endl;
U.write();
// Update/re-write phi
#include "createPhi.H"
phi.write();
singlePhaseTransportModel laminarTransport(U, phi);
autoPtr<incompressible::turbulenceModel> turbulence
(
incompressible::turbulenceModel::New(U, phi, laminarTransport)
);
if (isA<incompressible::RASModel>(turbulence()))
{
// Calculate nut - reference nut is calculated by the turbulence model
// on its construction
tmp<volScalarField> tnut = turbulence->nut();
volScalarField& nut = tnut();
volScalarField S(mag(dev(symm(fvc::grad(U)))));
nut = (1 - mask)*nut + mask*sqr(kappa*min(y, ybl))*::sqrt(2)*S;
// do not correct BC - wall functions will 'undo' manipulation above
// by using nut from turbulence model
if (args.optionFound("writenut"))
{
Info<< "Writing nut" << endl;
nut.write();
}
//--- Read and modify turbulence fields
// Turbulence k
tmp<volScalarField> tk = turbulence->k();
volScalarField& k = tk();
scalar ck0 = pow025(Cmu)*kappa;
k = (1 - mask)*k + mask*sqr(nut/(ck0*min(y, ybl)));
// do not correct BC - operation may use inconsistent fields wrt these
// local manipulations
// k.correctBoundaryConditions();
Info<< "Writing k\n" << endl;
k.write();
// Turbulence epsilon
tmp<volScalarField> tepsilon = turbulence->epsilon();
volScalarField& epsilon = tepsilon();
scalar ce0 = ::pow(Cmu, 0.75)/kappa;
epsilon = (1 - mask)*epsilon + mask*ce0*k*sqrt(k)/min(y, ybl);
// do not correct BC - wall functions will use non-updated k from
// turbulence model
// epsilon.correctBoundaryConditions();
Info<< "Writing epsilon\n" << endl;
epsilon.write();
// Turbulence omega
IOobject omegaHeader
(
"omega",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (omegaHeader.headerOk())
{
volScalarField omega(omegaHeader, mesh);
dimensionedScalar k0("VSMALL", k.dimensions(), VSMALL);
omega = (1 - mask)*omega + mask*epsilon/(Cmu*k + k0);
// do not correct BC - wall functions will use non-updated k from
// turbulence model
// omega.correctBoundaryConditions();
Info<< "Writing omega\n" << endl;
omega.write();
}
// Turbulence nuTilda
IOobject nuTildaHeader
(
"nuTilda",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
if (nuTildaHeader.headerOk())
{
volScalarField nuTilda(nuTildaHeader, mesh);
nuTilda = nut;
// do not correct BC
// nuTilda.correctBoundaryConditions();
Info<< "Writing nuTilda\n" << endl;
nuTilda.write();
}
}
Info<< nl << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //