Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / solvers / multiphase / settlingFoam / kEpsilon.H
blobac1de3b9e1c25152467a1f6ad04ae184e7f36bb4
1 if (turbulence)
3     if (mesh.changing())
4     {
5         y.correct();
6     }
8     dimensionedScalar kMin("kMin", k.dimensions(), SMALL);
9     dimensionedScalar epsilonMin("epsilonMin", epsilon.dimensions(), SMALL);
11     volScalarField divU(fvc::div(phi/fvc::interpolate(rho)));
13     tmp<volTensorField> tgradU = fvc::grad(U);
14     volScalarField G(2*mut*(tgradU() && dev(symm(tgradU()))));
15     tgradU.clear();
17     volScalarField Gcoef
18     (
19         Cmu*k/sigmak*(g & fvc::grad(rho))/(epsilon + epsilonMin)
20     );
22     #include "wallFunctions.H"
24     // Dissipation equation
25     fvScalarMatrix epsEqn
26     (
27         fvm::ddt(rho, epsilon)
28       + fvm::div(phi, epsilon)
29       - fvm::laplacian
30         (
31             mut/sigmaEps + mul, epsilon,
32             "laplacian(DepsilonEff,epsilon)"
33         )
34      ==
35         C1*G*epsilon/k
36       - fvm::SuSp(C1*(1.0 - C3)*Gcoef + (2.0/3.0*C1)*rho*divU, epsilon)
37       - fvm::Sp(C2*rho*epsilon/k, epsilon)
38     );
40     #include "wallDissipation.H"
42     epsEqn.relax();
43     epsEqn.solve();
45     bound(epsilon, epsilonMin);
48     // Turbulent kinetic energy equation
49     fvScalarMatrix kEqn
50     (
51         fvm::ddt(rho, k)
52       + fvm::div(phi, k)
53       - fvm::laplacian
54         (
55             mut/sigmak + mul, k,
56             "laplacian(DkEff,k)"
57         )
58      ==
59         G
60       - fvm::SuSp(Gcoef + 2.0/3.0*rho*divU, k)
61       - fvm::Sp(rho*epsilon/k, k)
62     );
64     kEqn.relax();
65     kEqn.solve();
67     bound(k, kMin);
70     //- Re-calculate viscosity
71     mut = rho*Cmu*sqr(k)/(epsilon + epsilonMin);
73     #include "wallViscosity.H"
76 mu = mut + mul;