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

210 lines
5.6 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2013 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"
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
template<class BasicTurbulenceModel>
void Foam::RASModel<BasicTurbulenceModel>::printCoeffs(const word& type)
{
if (printCoeffs_)
{
Info<< type << "Coeffs" << coeffDict_ << endl;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
Foam::RASModel<BasicTurbulenceModel>::RASModel
(
const word& type,
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
)
:
BasicTurbulenceModel
(
alpha,
rho,
U,
alphaPhi,
phi,
transport,
propertiesName
),
RASDict_(this->subOrEmptyDict("RAS")),
turbulence_(RASDict_.lookup("turbulence")),
printCoeffs_(RASDict_.lookupOrDefault<Switch>("printCoeffs", false)),
coeffDict_(RASDict_.subOrEmptyDict(type + "Coeffs")),
kMin_
(
dimensioned<scalar>::lookupOrAddToDict
(
"kMin",
RASDict_,
SMALL,
sqr(dimVelocity)
)
),
epsilonMin_
(
dimensioned<scalar>::lookupOrAddToDict
(
"epsilonMin",
RASDict_,
SMALL,
kMin_.dimensions()/dimTime
)
),
omegaMin_
(
dimensioned<scalar>::lookupOrAddToDict
(
"omegaMin",
RASDict_,
SMALL,
dimless/dimTime
)
)
{
// Force the construction of the mesh deltaCoeffs which may be needed
// for the construction of the derived models and BCs
this->mesh_.deltaCoeffs();
}
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
Foam::autoPtr<Foam::RASModel<BasicTurbulenceModel> >
Foam::RASModel<BasicTurbulenceModel>::New
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
)
{
// get model name, but do not register the dictionary
// otherwise it is registered in the database twice
const word modelType
(
IOdictionary
(
IOobject
(
IOobject::groupName(propertiesName, U.group()),
U.time().constant(),
U.db(),
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE,
false
)
).subDict("RAS").lookup("RASModel")
);
Info<< "Selecting RAS turbulence model " << modelType << endl;
typename dictionaryConstructorTable::iterator cstrIter =
dictionaryConstructorTablePtr_->find(modelType);
if (cstrIter == dictionaryConstructorTablePtr_->end())
{
FatalErrorIn
(
"RASModel::New"
"("
"const volScalarField&, "
"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()(alpha, rho, U, alphaPhi, phi, transport, propertiesName)
);
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
void Foam::RASModel<BasicTurbulenceModel>::correct()
{
BasicTurbulenceModel::correct();
}
template<class BasicTurbulenceModel>
bool Foam::RASModel<BasicTurbulenceModel>::read()
{
if (turbulenceModel::read())
{
RASDict_ <<= this->subDict("RAS");
RASDict_.lookup("turbulence") >> turbulence_;
if (const dictionary* dictPtr = RASDict_.subDictPtr(type() + "Coeffs"))
{
coeffDict_ <<= *dictPtr;
}
kMin_.readIfPresent(RASDict_);
epsilonMin_.readIfPresent(RASDict_);
omegaMin_.readIfPresent(RASDict_);
return true;
}
else
{
return false;
}
}
// ************************************************************************* //