Files
openfoam/applications/solvers/multiphase/bubbleFoam/createRASTurbulence.H
Mark Olesen 827e3659b9 consistency update: kappa=0.41, E=9.8
- this would be an argument for providing default values at the top-level
  compressible/incompressible turbulenceModel
2009-07-31 18:15:54 +02:00

179 lines
3.2 KiB
C

IOdictionary RASProperties
(
IOobject
(
"RASProperties",
runTime.constant(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
Switch turbulence
(
RASProperties.lookup("turbulence")
);
dictionary kEpsilonDict
(
RASProperties.subDictPtr("kEpsilonCoeffs")
);
dimensionedScalar Cmu
(
dimensionedScalar::lookupOrAddToDict
(
"Cmu",
kEpsilonDict,
0.09
)
);
dimensionedScalar C1
(
dimensionedScalar::lookupOrAddToDict
(
"C1",
kEpsilonDict,
1.44
)
);
dimensionedScalar C2
(
dimensionedScalar::lookupOrAddToDict
(
"C2",
kEpsilonDict,
1.92
)
);
dimensionedScalar alphak
(
dimensionedScalar::lookupOrAddToDict
(
"alphaEps",
kEpsilonDict,
1.0
)
);
dimensionedScalar alphaEps
(
dimensionedScalar::lookupOrAddToDict
(
"alphaEps",
kEpsilonDict,
0.76923
)
);
dictionary wallFunctionDict
(
RASProperties.subDictPtr("wallFunctionCoeffs")
);
dimensionedScalar kappa
(
dimensionedScalar::lookupOrAddToDict
(
"kappa",
wallFunctionDict,
0.41
)
);
dimensionedScalar E
(
dimensionedScalar::lookupOrAddToDict
(
"E",
wallFunctionDict,
9.8
)
);
if (RASProperties.lookupOrDefault<Switch>("printCoeffs", false))
{
Info<< "kEpsilonCoeffs" << kEpsilonDict << nl
<< "wallFunctionCoeffs" << wallFunctionDict << endl;
}
nearWallDist y(mesh);
Info<< "Reading field k\n" << endl;
volScalarField k
(
IOobject
(
"k",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
Info<< "Reading field epsilon\n" << endl;
volScalarField epsilon
(
IOobject
(
"epsilon",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
Info<< "Calculating field nutb\n" << endl;
volScalarField nutb
(
IOobject
(
"nutb",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
Cmu*sqr(k)/epsilon
);
Info<< "Calculating field nuEffa\n" << endl;
volScalarField nuEffa
(
IOobject
(
"nuEffa",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
sqr(Ct)*nutb + nua
);
Info<< "Calculating field nuEffb\n" << endl;
volScalarField nuEffb
(
IOobject
(
"nuEffb",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
nutb + nub
);