updating turb coeffs + fix typo

This commit is contained in:
andy
2009-08-13 09:54:57 +01:00
parent 62280a6e22
commit 9fb0ed08e7
2 changed files with 26 additions and 63 deletions

View File

@ -56,20 +56,11 @@
); );
dimensionedScalar rhoc dimensionedScalar rhoc(transportProperties.lookup("rhoc"));
(
transportProperties.lookup("rhoc")
);
dimensionedScalar rhod dimensionedScalar rhod(transportProperties.lookup("rhod"));
(
transportProperties.lookup("rhod")
);
dimensionedScalar muc dimensionedScalar muc(transportProperties.lookup("muc"));
(
transportProperties.lookup("muc")
);
dimensionedScalar plasticViscosityCoeff dimensionedScalar plasticViscosityCoeff
( (
@ -96,10 +87,7 @@
transportProperties.lookup("yieldStressOffset") transportProperties.lookup("yieldStressOffset")
); );
Switch BinghamPlastic Switch BinghamPlastic(transportProperties.lookup("BinghamPlastic"));
(
transportProperties.lookup("BinghamPlastic")
);
volScalarField rho volScalarField rho
( (
@ -147,8 +135,8 @@
IOobject::NO_READ, IOobject::NO_READ,
IOobject::AUTO_WRITE IOobject::AUTO_WRITE
), ),
muc + muc
plasticViscosity + plasticViscosity
( (
plasticViscosityCoeff, plasticViscosityCoeff,
plasticViscosityExponent, plasticViscosityExponent,
@ -176,10 +164,7 @@
Info<< "Selecting Drift-Flux model " << endl; Info<< "Selecting Drift-Flux model " << endl;
word VdjModel word VdjModel(transportProperties.lookup("VdjModel"));
(
transportProperties.lookup("VdjModel")
);
Info<< tab << VdjModel << " selected\n" << endl; Info<< tab << VdjModel << " selected\n" << endl;
@ -188,26 +173,13 @@
transportProperties.subDict(VdjModel + "Coeffs") transportProperties.subDict(VdjModel + "Coeffs")
); );
dimensionedVector V0 dimensionedVector V0(VdjModelCoeffs.lookup("V0"));
(
VdjModelCoeffs.lookup("V0")
);
dimensionedScalar a dimensionedScalar a(VdjModelCoeffs.lookup("a"));
(
VdjModelCoeffs.lookup("a")
);
dimensionedScalar a1 dimensionedScalar a1(VdjModelCoeffs.lookup("a1"));
(
VdjModelCoeffs.lookup("a1")
);
dimensionedScalar alphaMin
(
VdjModelCoeffs.lookup("alphaMin")
);
dimensionedScalar alphaMin(VdjModelCoeffs.lookup("alphaMin"));
IOdictionary RASProperties IOdictionary RASProperties
@ -223,15 +195,9 @@
); );
Switch turbulence Switch turbulence(RASProperties.lookup("turbulence"));
(
RASProperties.lookup("turbulence")
);
dictionary kEpsilonDict dictionary kEpsilonDict(RASProperties.subDictPtr("kEpsilonCoeffs"));
(
RASProperties.subDictPtr("kEpsilonCoeffs")
);
dimensionedScalar Cmu dimensionedScalar Cmu
( (
@ -273,30 +239,27 @@
) )
); );
dimensionedScalar alphak dimensionedScalar sigmak
( (
dimensionedScalar::lookupOrAddToDict dimensionedScalar::lookupOrAddToDict
( (
"alphaEps", "sigmak",
kEpsilonDict, kEpsilonDict,
1.0 1.0
) )
); );
dimensionedScalar alphaEps dimensionedScalar sigmaEps
( (
dimensionedScalar::lookupOrAddToDict dimensionedScalar::lookupOrAddToDict
( (
"alphaEps", "sigmaEps",
kEpsilonDict, kEpsilonDict,
0.76923 1.3
) )
); );
dictionary wallFunctionDict dictionary wallFunctionDict(RASProperties.subDictPtr("wallFunctionCoeffs"));
(
RASProperties.subDictPtr("wallFunctionCoeffs")
);
dimensionedScalar kappa dimensionedScalar kappa
( (

View File

@ -1,4 +1,4 @@
if(turbulence) if (turbulence)
{ {
if (mesh.changing()) if (mesh.changing())
{ {
@ -15,9 +15,9 @@ if(turbulence)
tgradU.clear(); tgradU.clear();
volScalarField Gcoef = volScalarField Gcoef =
alphak*Cmu*k*(g & fvc::grad(rho))/(epsilon + epsilon0); Cmu*k/sigmak*(g & fvc::grad(rho))/(epsilon + epsilon0);
# include "wallFunctions.H" #include "wallFunctions.H"
// Dissipation equation // Dissipation equation
fvScalarMatrix epsEqn fvScalarMatrix epsEqn
@ -26,7 +26,7 @@ if(turbulence)
+ fvm::div(phi, epsilon) + fvm::div(phi, epsilon)
- fvm::laplacian - fvm::laplacian
( (
alphaEps*mut + mul, epsilon, mut/sigmaEps + mul, epsilon,
"laplacian(DepsilonEff,epsilon)" "laplacian(DepsilonEff,epsilon)"
) )
== ==
@ -35,7 +35,7 @@ if(turbulence)
- fvm::Sp(C2*rho*epsilon/k, epsilon) - fvm::Sp(C2*rho*epsilon/k, epsilon)
); );
# include "wallDissipation.H" #include "wallDissipation.H"
epsEqn.relax(); epsEqn.relax();
epsEqn.solve(); epsEqn.solve();
@ -51,7 +51,7 @@ if(turbulence)
+ fvm::div(phi, k) + fvm::div(phi, k)
- fvm::laplacian - fvm::laplacian
( (
alphak*mut + mul, k, mut/sigmak + mul, k,
"laplacian(DkEff,k)" "laplacian(DkEff,k)"
) )
== ==
@ -66,7 +66,7 @@ if(turbulence)
//- Re-calculate viscosity //- Re-calculate viscosity
mut = rho*Cmu*sqr(k)/(epsilon + epsilon0); mut = rho*Cmu*sqr(k)/(epsilon + epsilon0);
# include "wallViscosity.H" #include "wallViscosity.H"
} }
mu = mut + mul; mu = mut + mul;