mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
198 lines
5.2 KiB
C
198 lines
5.2 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2013-2017 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<< coeffDict_.dictName() << coeffDict_ << endl;
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
template<class BasicTurbulenceModel>
|
|
Foam::RASModel<BasicTurbulenceModel>::RASModel
|
|
(
|
|
const word& type,
|
|
const alphaField& alpha,
|
|
const rhoField& rho,
|
|
const volVectorField& U,
|
|
const surfaceScalarField& alphaRhoPhi,
|
|
const surfaceScalarField& phi,
|
|
const transportModel& transport,
|
|
const word& propertiesName
|
|
)
|
|
:
|
|
BasicTurbulenceModel
|
|
(
|
|
type,
|
|
alpha,
|
|
rho,
|
|
U,
|
|
alphaRhoPhi,
|
|
phi,
|
|
transport,
|
|
propertiesName
|
|
),
|
|
|
|
RASDict_(this->subOrEmptyDict("RAS")),
|
|
turbulence_(RASDict_.get<Switch>("turbulence")),
|
|
printCoeffs_(RASDict_.lookupOrDefault<Switch>("printCoeffs", false)),
|
|
coeffDict_(RASDict_.optionalSubDict(type + "Coeffs")),
|
|
|
|
kMin_
|
|
(
|
|
dimensioned<scalar>::lookupOrAddToDict
|
|
(
|
|
"kMin",
|
|
RASDict_,
|
|
sqr(dimVelocity),
|
|
SMALL
|
|
)
|
|
),
|
|
|
|
epsilonMin_
|
|
(
|
|
dimensioned<scalar>::lookupOrAddToDict
|
|
(
|
|
"epsilonMin",
|
|
RASDict_,
|
|
kMin_.dimensions()/dimTime,
|
|
SMALL
|
|
)
|
|
),
|
|
|
|
omegaMin_
|
|
(
|
|
dimensioned<scalar>::lookupOrAddToDict
|
|
(
|
|
"omegaMin",
|
|
RASDict_,
|
|
dimless/dimTime,
|
|
SMALL
|
|
)
|
|
)
|
|
{
|
|
// 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& alphaRhoPhi,
|
|
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, alphaRhoPhi.group()),
|
|
U.time().constant(),
|
|
U.db(),
|
|
IOobject::MUST_READ_IF_MODIFIED,
|
|
IOobject::NO_WRITE,
|
|
false
|
|
)
|
|
).subDict("RAS").get<word>("RASModel")
|
|
);
|
|
|
|
Info<< "Selecting RAS turbulence model " << modelType << endl;
|
|
|
|
auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
|
|
|
|
if (!cstrIter.found())
|
|
{
|
|
FatalErrorInFunction
|
|
<< "Unknown RASModel type "
|
|
<< modelType << nl << nl
|
|
<< "Valid RASModel types:" << endl
|
|
<< dictionaryConstructorTablePtr_->sortedToc()
|
|
<< exit(FatalError);
|
|
}
|
|
|
|
return autoPtr<RASModel>
|
|
(
|
|
cstrIter()(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
|
|
);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
template<class BasicTurbulenceModel>
|
|
bool Foam::RASModel<BasicTurbulenceModel>::read()
|
|
{
|
|
if (BasicTurbulenceModel::read())
|
|
{
|
|
RASDict_ <<= this->subDict("RAS");
|
|
RASDict_.readEntry("turbulence", turbulence_);
|
|
|
|
coeffDict_ <<= RASDict_.optionalSubDict(type() + "Coeffs");
|
|
|
|
kMin_.readIfPresent(RASDict_);
|
|
epsilonMin_.readIfPresent(RASDict_);
|
|
omegaMin_.readIfPresent(RASDict_);
|
|
|
|
return true;
|
|
}
|
|
else
|
|
{
|
|
return false;
|
|
}
|
|
}
|
|
|
|
|
|
template<class BasicTurbulenceModel>
|
|
void Foam::RASModel<BasicTurbulenceModel>::correct()
|
|
{
|
|
BasicTurbulenceModel::correct();
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|