Files
openfoam/src/turbulenceModels/RAS/incompressible/kEpsilon/kEpsilon.C
2008-09-29 10:14:41 +01:00

282 lines
6.1 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2008 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 "kEpsilon.H"
#include "addToRunTimeSelectionTable.H"
#include "wallFvPatch.H"
#include "backwardsCompatibilityWallFunctions.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(kEpsilon, 0);
addToRunTimeSelectionTable(RASModel, kEpsilon, dictionary);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
// Construct from components
kEpsilon::kEpsilon
(
const volVectorField& U,
const surfaceScalarField& phi,
transportModel& lamTransportModel
)
:
RASModel(typeName, U, phi, lamTransportModel),
Cmu_
(
dimensioned<scalar>::lookupOrAddToDict
(
"Cmu",
coeffDict_,
0.09
)
),
C1_
(
dimensioned<scalar>::lookupOrAddToDict
(
"C1",
coeffDict_,
1.44
)
),
C2_
(
dimensioned<scalar>::lookupOrAddToDict
(
"C2",
coeffDict_,
1.92
)
),
alphaEps_
(
dimensioned<scalar>::lookupOrAddToDict
(
"alphaEps",
coeffDict_,
0.76923
)
),
k_
(
IOobject
(
"k",
runTime_.timeName(),
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
autoCreateK("k", mesh_)
),
epsilon_
(
IOobject
(
"epsilon",
runTime_.timeName(),
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
autoCreateEpsilon("epsilon", mesh_)
),
nut_
(
IOobject
(
"nut",
runTime_.timeName(),
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
autoCreateNut("nut", mesh_)
)
{
nut_ == Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_);
nut_.correctBoundaryConditions();
printCoeffs();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
tmp<volSymmTensorField> kEpsilon::R() const
{
return tmp<volSymmTensorField>
(
new volSymmTensorField
(
IOobject
(
"R",
runTime_.timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
k_.boundaryField().types()
)
);
}
tmp<volSymmTensorField> kEpsilon::devReff() const
{
return tmp<volSymmTensorField>
(
new volSymmTensorField
(
IOobject
(
"devRhoReff",
runTime_.timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
-nuEff()*dev(twoSymm(fvc::grad(U_)))
)
);
}
tmp<fvVectorMatrix> kEpsilon::divDevReff(volVectorField& U) const
{
return
(
- fvm::laplacian(nuEff(), U)
- fvc::div(nuEff()*dev(fvc::grad(U)().T()))
);
}
bool kEpsilon::read()
{
if (RASModel::read())
{
Cmu_.readIfPresent(coeffDict_);
C1_.readIfPresent(coeffDict_);
C2_.readIfPresent(coeffDict_);
alphaEps_.readIfPresent(coeffDict_);
return true;
}
else
{
return false;
}
}
void kEpsilon::correct()
{
transportModel_.correct();
if (!turbulence_)
{
return;
}
RASModel::correct();
if (mesh_.changing())
{
y_.correct();
}
volScalarField G("G", nut_*2*magSqr(symm(fvc::grad(U_))));
// Update espsilon and G at the wall
epsilon_.boundaryField().updateCoeffs();
// Dissipation equation
tmp<fvScalarMatrix> epsEqn
(
fvm::ddt(epsilon_)
+ fvm::div(phi_, epsilon_)
- fvm::Sp(fvc::div(phi_), epsilon_)
- fvm::laplacian(DepsilonEff(), epsilon_)
==
C1_*G*epsilon_/k_
- fvm::Sp(C2_*epsilon_/k_, epsilon_)
);
epsEqn().relax();
epsEqn().boundaryManipulate(epsilon_.boundaryField());
solve(epsEqn);
bound(epsilon_, epsilon0_);
// Turbulent kinetic energy equation
tmp<fvScalarMatrix> kEqn
(
fvm::ddt(k_)
+ fvm::div(phi_, k_)
- fvm::Sp(fvc::div(phi_), k_)
- fvm::laplacian(DkEff(), k_)
==
G
- fvm::Sp(epsilon_/k_, k_)
);
kEqn().relax();
solve(kEqn);
bound(k_, k0_);
// Re-calculate viscosity
nut_ == Cmu_*sqr(k_)/epsilon_;
nut_.correctBoundaryConditions();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// ************************************************************************* //