mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
207 lines
5.5 KiB
C
207 lines
5.5 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/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "RASModel.H"
|
|
#include "addToRunTimeSelectionTable.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
namespace incompressible
|
|
{
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
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 volVectorField& U,
|
|
const surfaceScalarField& phi,
|
|
const transportModel& transport,
|
|
const word& turbulenceModelName
|
|
)
|
|
:
|
|
turbulenceModel(U, phi, transport, turbulenceModelName),
|
|
|
|
IOdictionary
|
|
(
|
|
IOobject
|
|
(
|
|
"RASProperties",
|
|
U.time().constant(),
|
|
U.db(),
|
|
IOobject::MUST_READ_IF_MODIFIED,
|
|
IOobject::NO_WRITE
|
|
)
|
|
),
|
|
|
|
turbulence_(lookup("turbulence")),
|
|
printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
|
|
coeffDict_(subOrEmptyDict(type + "Coeffs")),
|
|
|
|
kMin_("kMin", sqr(dimVelocity), SMALL),
|
|
epsilonMin_("epsilonMin", kMin_.dimensions()/dimTime, SMALL),
|
|
omegaMin_("omegaMin", dimless/dimTime, SMALL)
|
|
{
|
|
kMin_.readIfPresent(*this);
|
|
epsilonMin_.readIfPresent(*this);
|
|
omegaMin_.readIfPresent(*this);
|
|
|
|
// 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 volVectorField& U,
|
|
const surfaceScalarField& phi,
|
|
const transportModel& transport,
|
|
const word& turbulenceModelName
|
|
)
|
|
{
|
|
// get model name, but do not register the dictionary
|
|
// otherwise it is registered in the database twice
|
|
const word modelType
|
|
(
|
|
IOdictionary
|
|
(
|
|
IOobject
|
|
(
|
|
"RASProperties",
|
|
U.time().constant(),
|
|
U.db(),
|
|
IOobject::MUST_READ_IF_MODIFIED,
|
|
IOobject::NO_WRITE,
|
|
false
|
|
)
|
|
).lookup("RASModel")
|
|
);
|
|
|
|
Info<< "Selecting RAS turbulence model " << modelType << endl;
|
|
|
|
dictionaryConstructorTable::iterator cstrIter =
|
|
dictionaryConstructorTablePtr_->find(modelType);
|
|
|
|
if (cstrIter == dictionaryConstructorTablePtr_->end())
|
|
{
|
|
FatalErrorIn
|
|
(
|
|
"RASModel::New"
|
|
"("
|
|
"const volVectorField&, "
|
|
"const surfaceScalarField&, "
|
|
"transportModel&, "
|
|
"const word&"
|
|
")"
|
|
) << "Unknown RASModel type "
|
|
<< modelType << nl << nl
|
|
<< "Valid RASModel types:" << endl
|
|
<< dictionaryConstructorTablePtr_->sortedToc()
|
|
<< exit(FatalError);
|
|
}
|
|
|
|
return autoPtr<RASModel>
|
|
(
|
|
cstrIter()(U, phi, transport, turbulenceModelName)
|
|
);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void RASModel::correct()
|
|
{
|
|
turbulenceModel::correct();
|
|
}
|
|
|
|
|
|
bool RASModel::read()
|
|
{
|
|
//if (regIOobject::read())
|
|
|
|
// Bit of trickery : we are both IOdictionary ('RASProperties') and
|
|
// an regIOobject from the turbulenceModel level. Problem is to distinguish
|
|
// between the two - we only want to reread the IOdictionary.
|
|
|
|
bool ok = IOdictionary::readData
|
|
(
|
|
IOdictionary::readStream
|
|
(
|
|
IOdictionary::type()
|
|
)
|
|
);
|
|
IOdictionary::close();
|
|
|
|
if (ok)
|
|
{
|
|
lookup("turbulence") >> turbulence_;
|
|
|
|
if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
|
|
{
|
|
coeffDict_ <<= *dictPtr;
|
|
}
|
|
|
|
kMin_.readIfPresent(*this);
|
|
epsilonMin_.readIfPresent(*this);
|
|
omegaMin_.readIfPresent(*this);
|
|
|
|
return true;
|
|
}
|
|
else
|
|
{
|
|
return false;
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
} // End namespace incompressible
|
|
} // End namespace Foam
|
|
|
|
// ************************************************************************* //
|