Files
openfoam/src/TurbulenceModels/turbulenceModels/laminar/laminarModel/laminarModel.C
Kutalmis Bercin d7622d1ef8 ENH: add virtual omega to turbulenceModels
STYLE: update dimensionSets in SpalartAllmaras
2020-03-11 16:31:42 +00:00

336 lines
8.4 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2016-2017 OpenFOAM Foundation
Copyright (C) 2019-2020 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "laminarModel.H"
#include "Stokes.H"
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
template<class BasicTurbulenceModel>
void Foam::laminarModel<BasicTurbulenceModel>::printCoeffs(const word& type)
{
if (printCoeffs_)
{
Info<< coeffDict_.dictName() << coeffDict_ << endl;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
Foam::laminarModel<BasicTurbulenceModel>::laminarModel
(
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
),
laminarDict_(this->subOrEmptyDict("laminar")),
printCoeffs_(laminarDict_.lookupOrDefault<Switch>("printCoeffs", false)),
coeffDict_(laminarDict_.optionalSubDict(type + "Coeffs"))
{
// 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::laminarModel<BasicTurbulenceModel>>
Foam::laminarModel<BasicTurbulenceModel>::New
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName
)
{
const IOdictionary modelDict
(
IOobject
(
IOobject::groupName(propertiesName, alphaRhoPhi.group()),
U.time().constant(),
U.db(),
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE,
false // Do not register
)
);
const dictionary* dictptr = modelDict.findDict("laminar");
if (dictptr)
{
const dictionary& dict = *dictptr;
const word modelType(dict.get<word>("laminarModel"));
Info<< "Selecting laminar stress model " << modelType << endl;
auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
if (!cstrIter.found())
{
FatalIOErrorInLookup
(
dict,
"laminarModel",
modelType,
*dictionaryConstructorTablePtr_
) << exit(FatalIOError);
}
return autoPtr<laminarModel>
(
cstrIter()
(
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport, propertiesName)
);
}
else
{
Info<< "Selecting laminar stress model "
<< laminarModels::Stokes<BasicTurbulenceModel>::typeName << endl;
return autoPtr<laminarModel>
(
new laminarModels::Stokes<BasicTurbulenceModel>
(
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport,
propertiesName
)
);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class BasicTurbulenceModel>
bool Foam::laminarModel<BasicTurbulenceModel>::read()
{
if (BasicTurbulenceModel::read())
{
laminarDict_ <<= this->subDict("laminar");
coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
return true;
}
return false;
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volScalarField>
Foam::laminarModel<BasicTurbulenceModel>::nut() const
{
return tmp<volScalarField>::New
(
IOobject
(
IOobject::groupName("nut", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh_,
dimensionedScalar(dimViscosity, Zero)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::scalarField>
Foam::laminarModel<BasicTurbulenceModel>::nut
(
const label patchi
) const
{
return tmp<scalarField>
(
new scalarField(this->mesh_.boundary()[patchi].size(), Zero)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volScalarField>
Foam::laminarModel<BasicTurbulenceModel>::nuEff() const
{
return tmp<volScalarField>
(
new volScalarField
(
IOobject::groupName("nuEff", this->alphaRhoPhi_.group()), this->nu()
)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::scalarField>
Foam::laminarModel<BasicTurbulenceModel>::nuEff
(
const label patchi
) const
{
return this->nu(patchi);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volScalarField>
Foam::laminarModel<BasicTurbulenceModel>::k() const
{
return tmp<volScalarField>::New
(
IOobject
(
IOobject::groupName("k", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh_,
dimensionedScalar(sqr(this->U_.dimensions()), Zero)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volScalarField>
Foam::laminarModel<BasicTurbulenceModel>::epsilon() const
{
return tmp<volScalarField>::New
(
IOobject
(
IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh_,
dimensionedScalar(sqr(this->U_.dimensions())/dimTime, Zero)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volScalarField>
Foam::laminarModel<BasicTurbulenceModel>::omega() const
{
return tmp<volScalarField>::New
(
IOobject
(
IOobject::groupName("omega", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh_,
dimensionedScalar(dimless/dimTime, Zero)
);
}
template<class BasicTurbulenceModel>
Foam::tmp<Foam::volSymmTensorField>
Foam::laminarModel<BasicTurbulenceModel>::R() const
{
return tmp<volSymmTensorField>::New
(
IOobject
(
IOobject::groupName("R", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
this->mesh_,
dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
);
}
template<class BasicTurbulenceModel>
void Foam::laminarModel<BasicTurbulenceModel>::correct()
{
BasicTurbulenceModel::correct();
}
// ************************************************************************* //