Files
openfoam/src/turbulenceModels/compressible/RAS/RASModel/RASModel.C

243 lines
6.3 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "RASModel.H"
#include "wallFvPatch.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace compressible
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(RASModel, 0);
defineRunTimeSelectionTable(RASModel, dictionary);
addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
void RASModel::printCoeffs()
{
if (printCoeffs_)
{
Info<< type() << "Coeffs" << coeffDict_ << endl;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
RASModel::RASModel
(
const word& type,
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const basicThermo& thermophysicalModel
)
:
turbulenceModel(rho, U, phi, thermophysicalModel),
IOdictionary
(
IOobject
(
"RASProperties",
U.time().constant(),
U.db(),
IOobject::MUST_READ,
IOobject::NO_WRITE
)
),
turbulence_(lookup("turbulence")),
printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
coeffDict_(subDictPtr(type + "Coeffs")),
k0_("k0", dimVelocity*dimVelocity, SMALL),
epsilon0_("epsilon", k0_.dimensions()/dimTime, SMALL),
epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
omega0_("omega", dimless/dimTime, SMALL),
omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
y_(mesh_)
{
// Force the construction of the mesh deltaCoeffs which may be needed
// for the construction of the derived models and BCs
mesh_.deltaCoeffs();
}
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
autoPtr<RASModel> RASModel::New
(
const volScalarField& rho,
const volVectorField& U,
const surfaceScalarField& phi,
const basicThermo& thermophysicalModel
)
{
word modelName;
// Enclose the creation of the dictionary to ensure it is deleted
// before the turbulenceModel is created otherwise the dictionary is
// entered in the database twice
{
IOdictionary dict
(
IOobject
(
"RASProperties",
U.time().constant(),
U.db(),
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
dict.lookup("RASModel") >> modelName;
}
Info<< "Selecting RAS turbulence model " << modelName << endl;
dictionaryConstructorTable::iterator cstrIter =
dictionaryConstructorTablePtr_->find(modelName);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorIn
(
"RASModel::New(const volScalarField&, "
"const volVectorField&, const surfaceScalarField&, "
"basicThermo&)"
) << "Unknown RASModel type " << modelName
<< endl << endl
<< "Valid RASModel types are :" << endl
<< dictionaryConstructorTablePtr_->sortedToc()
<< exit(FatalError);
}
return autoPtr<RASModel>
(
cstrIter()(rho, U, phi, thermophysicalModel)
);
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
{
scalar ypl = 11.0;
for (int i=0; i<10; i++)
{
ypl = log(E*ypl)/kappa;
}
return ypl;
}
tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
{
const fvPatch& curPatch = mesh_.boundary()[patchNo];
tmp<scalarField> tYp(new scalarField(curPatch.size()));
scalarField& Yp = tYp();
if (isA<wallFvPatch>(curPatch))
{
Yp = pow025(Cmu)
*y_[patchNo]
*sqrt(k()().boundaryField()[patchNo].patchInternalField())
/(
mu().boundaryField()[patchNo].patchInternalField()
/rho_.boundaryField()[patchNo]
);
}
else
{
WarningIn
(
"tmp<scalarField> RASModel::yPlus(const label patchNo) const"
) << "Patch " << patchNo << " is not a wall. Returning null field"
<< nl << endl;
Yp.setSize(0);
}
return tYp;
}
void RASModel::correct()
{
if (mesh_.changing())
{
y_.correct();
}
}
bool RASModel::read()
{
if (regIOobject::read())
{
lookup("turbulence") >> turbulence_;
if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
{
coeffDict_ <<= *dictPtr;
}
k0_.readIfPresent(*this);
epsilon0_.readIfPresent(*this);
epsilonSmall_.readIfPresent(*this);
omega0_.readIfPresent(*this);
omegaSmall_.readIfPresent(*this);
return true;
}
else
{
return false;
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace compressible
} // End namespace Foam
// ************************************************************************* //