tmp<volScalarField> realizableKE::rCmu ( const volTensorField& gradU ) { volScalarField S2 = 2*magSqr(dev(symm(gradU))); volScalarField magS = sqrt(S2); return rCmu(gradU, S2, magS); }
realizableKE::realizableKE ( const volVectorField& U, const surfaceScalarField& phi, transportModel& lamTransportModel ) : RASModel(typeName, U, phi, lamTransportModel), Cmu_ ( dimensioned<scalar>::lookupOrAddToDict ( "Cmu", coeffDict_, 0.09 ) ), A0_ ( dimensioned<scalar>::lookupOrAddToDict ( "A0", coeffDict_, 4.0 ) ), C2_ ( dimensioned<scalar>::lookupOrAddToDict ( "C2", coeffDict_, 1.9 ) ), sigmak_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmak", coeffDict_, 1.0 ) ), sigmaEps_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmaEps", coeffDict_, 1.2 ) ), 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_) ) { bound(k_, k0_); bound(epsilon_, epsilon0_); nut_ = rCmu(fvc::grad(U_))*sqr(k_)/(epsilon_ + epsilonSmall_); nut_.correctBoundaryConditions(); printCoeffs(); }
realizableKE_Veh::realizableKE_Veh ( const volVectorField& U, const surfaceScalarField& phi, transportModel& transport, const word& turbulenceModelName, const word& modelName ) : RASModel(modelName, U, phi, transport, turbulenceModelName), Cfcar_ ( dimensioned<scalar>::lookupOrAddToDict ( "Cfcar", coeffDict_, 0.64 ) ), L_ ( dimensioned<scalar>::lookupOrAddToDict ( "L", coeffDict_, 10, dimensionSet(0,1,0,0,0,0,0) ) ), Cecar_ ( dimensioned<scalar>::lookupOrAddToDict ( "Cecar", coeffDict_, 2.47 ) ), Cmu_ ( dimensioned<scalar>::lookupOrAddToDict ( "Cmu", coeffDict_, 0.09 ) ), A0_ ( dimensioned<scalar>::lookupOrAddToDict ( "A0", coeffDict_, 4.0 ) ), C2_ ( dimensioned<scalar>::lookupOrAddToDict ( "C2", coeffDict_, 1.9 ) ), sigmak_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmak", coeffDict_, 1.0 ) ), sigmaEps_ ( dimensioned<scalar>::lookupOrAddToDict ( "sigmaEps", coeffDict_, 1.2 ) ), 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_) ), Ucar_ ( IOobject ( "Ucar", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ), VAD_ ( IOobject ( "VAD", runTime_.timeName(), mesh_, IOobject::MUST_READ, IOobject::AUTO_WRITE ), mesh_ ) { bound(k_, kMin_); bound(epsilon_, epsilonMin_); nut_ = rCmu(fvc::grad(U_))*sqr(k_)/epsilon_; nut_.correctBoundaryConditions(); printCoeffs(); }
void realizableKE::correct() { // Bound in case of topological change // HJ, 22/Aug/2007 if (mesh_.changing()) { bound(k_, k0_); bound(epsilon_, epsilon0_); } RASModel::correct(); if (!turbulence_) { return; } volTensorField gradU = fvc::grad(U_); volScalarField S2 = 2*magSqr(dev(symm(gradU))); volScalarField magS = sqrt(S2); volScalarField eta = magS*k_/epsilon_; volScalarField C1 = max(eta/(scalar(5) + eta), scalar(0.43)); volScalarField G("RASModel::G", nut_*S2); // Update epsilon and G at the wall epsilon_.boundaryField().updateCoeffs(); // Dissipation equation tmp<fvScalarMatrix> epsEqn ( fvm::ddt(epsilon_) + fvm::div(phi_, epsilon_) + fvm::SuSp(-fvc::div(phi_), epsilon_) - fvm::laplacian(DepsilonEff(), epsilon_) == C1*magS*epsilon_ - fvm::Sp ( C2_*epsilon_/(k_ + sqrt(nu()*epsilon_)), 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::SuSp(-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_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_; nut_.correctBoundaryConditions(); }
void realizableKE_Veh::correct() { RASModel::correct(); if (!turbulence_) { return; } const volTensorField gradU(fvc::grad(U_)); const volScalarField S2(2*magSqr(dev(symm(gradU)))); const volScalarField magS(sqrt(S2)); const volScalarField eta(magS*k_/epsilon_); tmp<volScalarField> C1 = max(eta/(scalar(5) + eta), scalar(0.43)); volScalarField G("RASModel::G", nut_*S2); //Estimating source terms for vehicle canopy volVectorField Fi = 0.5*Cfcar_*VAD_*mag(U_-Ucar_)*(U_-Ucar_); volScalarField Fk = (U_-Ucar_)&Fi; volScalarField Fepsilon = epsilon_*sqrt(k_)*Cecar_/L_; // Update epsilon 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_) == Fepsilon + C1*magS*epsilon_ - fvm::Sp ( C2_*epsilon_/(k_ + sqrt(nu()*epsilon_)), epsilon_ ) ); epsEqn().relax(); epsEqn().boundaryManipulate(epsilon_.boundaryField()); solve(epsEqn); bound(epsilon_, epsilonMin_); // Turbulent kinetic energy equation tmp<fvScalarMatrix> kEqn ( fvm::ddt(k_) + fvm::div(phi_, k_) - fvm::Sp(fvc::div(phi_), k_) - fvm::laplacian(DkEff(), k_) == Fk + G - fvm::Sp(epsilon_/k_, k_) ); kEqn().relax(); solve(kEqn); bound(k_, kMin_); // Re-calculate viscosity nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_; nut_.correctBoundaryConditions(); }