LaunderSharmaKE::LaunderSharmaKE ( 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 ) ), sigmaEps_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmaEps", coeffDict_, 1.3 ) ), k_ ( IOobject ( "k", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), epsilonTilda_ ( IOobject ( "epsilon", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), nut_ ( IOobject ( "nut", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateLowReNut("nut", mesh_) ) { nut_ = Cmu_*fMu()*sqr(k_)/(epsilonTilda_ + epsilonSmall_); nut_.correctBoundaryConditions(); printCoeffs(); }
void LaunderSharmaKE::correct() { // Bound in case of topological change // HJ, 22/Aug/2007 if (mesh_.changing()) { bound(k_, k0_); bound(epsilonTilda_, epsilon0_); } RASModel::correct(); if (!turbulence_) { return; } volScalarField S2 = 2*magSqr(symm(fvc::grad(U_))); volScalarField G("RASModel::G", nut_*S2); volScalarField E = 2.0*nu()*nut_*fvc::magSqrGradGrad(U_); volScalarField D = 2.0*nu()*magSqr(fvc::grad(sqrt(k_))); // Dissipation rate equation tmp<fvScalarMatrix> epsEqn ( fvm::ddt(epsilonTilda_) + fvm::div(phi_, epsilonTilda_) + fvm::SuSp(-fvc::div(phi_), epsilonTilda_) - fvm::laplacian(DepsilonEff(), epsilonTilda_) == C1_*G*epsilonTilda_/k_ - fvm::Sp(C2_*f2()*epsilonTilda_/k_, epsilonTilda_) + E ); epsEqn().relax(); solve(epsEqn); bound(epsilonTilda_, epsilon0_); // Turbulent kinetic energy equation tmp<fvScalarMatrix> kEqn ( fvm::ddt(k_) + fvm::div(phi_, k_) + fvm::SuSp(-fvc::div(phi_), k_) - fvm::laplacian(DkEff(), k_) == G - fvm::Sp((epsilonTilda_ + D)/k_, k_) ); kEqn().relax(); solve(kEqn); bound(k_, k0_); // Re-calculate viscosity nut_ == Cmu_*fMu()*sqr(k_)/epsilonTilda_; }
LaunderSharmaKE::LaunderSharmaKE ( const volScalarField& rho, const volVectorField& U, const surfaceScalarField& phi, const basicThermo& thermophysicalModel ) : RASModel(typeName, rho, U, phi, thermophysicalModel), Cmu_ ( dimensionedScalar::lookupOrAddToDict ( "Cmu", coeffDict_, 0.09 ) ), C1_ ( dimensionedScalar::lookupOrAddToDict ( "C1", coeffDict_, 1.44 ) ), C2_ ( dimensionedScalar::lookupOrAddToDict ( "C2", coeffDict_, 1.92 ) ), C3_ ( dimensionedScalar::lookupOrAddToDict ( "C3", coeffDict_, -0.33 ) ), sigmak_ ( dimensionedScalar::lookupOrAddToDict ( "sigmak", coeffDict_, 1.0 ) ), sigmaEps_ ( dimensionedScalar::lookupOrAddToDict ( "sigmaEps", coeffDict_, 1.3 ) ), Prt_ ( dimensionedScalar::lookupOrAddToDict ( "Prt", coeffDict_, 1.0 ) ), k_ ( IOobject ( "k", runTime_.timeName(), U_.db(), IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), epsilon_ ( IOobject ( "epsilon", runTime_.timeName(), U_.db(), IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), mut_ ( IOobject ( "mut", runTime_.timeName(), U_.db(), IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateLowReMut("mut", mesh_) ), alphat_ ( IOobject ( "alphat", runTime_.timeName(), U_.db(), IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateAlphat("alphat", mesh_) ) { mut_ = rho_*Cmu_*fMu()*sqr(k_)/(epsilon_ + epsilonSmall_); mut_ = min(mut_, muRatio()*mu()); mut_.correctBoundaryConditions(); alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); printCoeffs(); }
void LaunderSharmaKE::correct() { // Bound in case of topological change // HJ, 22/Aug/2007 if (mesh_.changing()) { bound(k_, k0_); bound(epsilon_, epsilon0_); } if (!turbulence_) { // Re-calculate viscosity mut_ == rho_*Cmu_*fMu()*sqr(k_)/(epsilon_ + epsilonSmall_); mut_ == min(mut_, muRatio()*mu()); // Re-calculate thermal diffusivity alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); return; } RASModel::correct(); // Calculate parameters and coefficients for Launder-Sharma low-Reynolds // number model volScalarField E = 2.0*mu()*mut_*fvc::magSqrGradGrad(U_)/rho_; volScalarField D = 2.0*mu()*magSqr(fvc::grad(sqrt(k_)))/rho_; volScalarField divU = fvc::div(phi_/fvc::interpolate(rho_)); if (mesh_.moving()) { divU += fvc::div(mesh_.phi()); } tmp<volTensorField> tgradU = fvc::grad(U_); volScalarField G("RASModel::G", mut_*(tgradU() && dev(twoSymm(tgradU())))); tgradU.clear(); // Dissipation equation tmp<fvScalarMatrix> epsEqn ( fvm::ddt(rho_, epsilon_) + fvm::div(phi_, epsilon_) - fvm::laplacian(DepsilonEff(), epsilon_) == C1_*G*epsilon_/k_ + fvm::SuSp((C3_ - 2.0/3.0*C1_)*rho_*divU, epsilon_) - fvm::Sp(C2_*f2()*rho_*epsilon_/k_, epsilon_) //+ 0.75*1.5*flameKproduction*epsilon_/k_ + E ); epsEqn().relax(); solve(epsEqn); bound(epsilon_, epsilon0_); // Turbulent kinetic energy equation tmp<fvScalarMatrix> kEqn ( fvm::ddt(rho_, k_) + fvm::div(phi_, k_) - fvm::laplacian(DkEff(), k_) == G - fvm::SuSp(2.0/3.0*rho_*divU, k_) - fvm::Sp(rho_*(epsilon_ + D)/k_, k_) //+ flameKproduction ); kEqn().relax(); solve(kEqn); bound(k_, k0_); // Re-calculate viscosity mut_ == Cmu_*fMu()*rho_*sqr(k_)/(epsilon_ + epsilonSmall_); mut_ == min(mut_, muRatio()*mu()); // Re-calculate thermal diffusivity alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); }
qZeta::qZeta ( 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 ) ), sigmaZeta_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmaZeta", coeffDict_, 1.3 ) ), anisotropic_ ( Switch::lookupOrAddToDict ( "anisotropic", coeffDict_, false ) ), k_ ( IOobject ( "k", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), epsilon_ ( IOobject ( "epsilon", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), q_ ( IOobject ( "q", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), sqrt(k_), k_.boundaryField().types() ), zeta_ ( IOobject ( "zeta", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), epsilon_/(2.0*q_), epsilon_.boundaryField().types() ), nut_ ( IOobject ( "nut", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateNut("nut", mesh_) ) { nut_ = Cmu_*fMu()*sqr(k_)/(epsilon_ + epsilonSmall_); nut_.correctBoundaryConditions(); printCoeffs(); }
void qZeta::correct() { RASModel::correct(); if (!turbulence_) { return; } volScalarField S2 = 2*magSqr(symm(fvc::grad(U_))); volScalarField G("RASModel::G", nut_/(2.0*q_)*S2); volScalarField E = nu()*nut_/q_*fvc::magSqrGradGrad(U_); // Zeta equation tmp<fvScalarMatrix> zetaEqn ( fvm::ddt(zeta_) + fvm::div(phi_, zeta_) - fvm::laplacian(DzetaEff(), zeta_) == (2.0*C1_ - 1)*G*zeta_/q_ - fvm::Sp((2.0*C2_ - dimensionedScalar(1.0))*f2()*zeta_/q_, zeta_) + E ); zetaEqn().relax(); solve(zetaEqn); bound(zeta_, epsilon0_/(2*sqrt(k0_))); // q equation tmp<fvScalarMatrix> qEqn ( fvm::ddt(q_) + fvm::div(phi_, q_) - fvm::laplacian(DqEff(), q_) == G - fvm::Sp(zeta_/q_, q_) ); qEqn().relax(); solve(qEqn); bound(q_, sqrt(k0_)); // Re-calculate k and epsilon k_ = sqr(q_); k_.correctBoundaryConditions(); epsilon_ = 2*q_*zeta_; epsilon_.correctBoundaryConditions(); // Re-calculate viscosity nut_ = Cmu_*fMu()*sqr(k_)/epsilon_; nut_.correctBoundaryConditions(); }