/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 3.2 \\ / A nd | Web: http://www.foam-extend.org \\/ M anipulation | For copyright notice see file Copyright ------------------------------------------------------------------------------- License This file is part of foam-extend. foam-extend is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. foam-extend is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with foam-extend. If not, see . \*---------------------------------------------------------------------------*/ #include "PDRkEpsilon.H" #include "PDRDragModel.H" #include "addToRunTimeSelectionTable.H" #include "backwardsCompatibilityWallFunctions.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { namespace compressible { namespace RASModels { // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // defineTypeNameAndDebug(PDRkEpsilon, 0); addToRunTimeSelectionTable(RASModel, PDRkEpsilon, dictionary); // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // PDRkEpsilon::PDRkEpsilon ( const volScalarField& rho, const volVectorField& U, const surfaceScalarField& phi, const basicThermo& thermophysicalModel ) : RASModel(typeName, rho, U, phi, thermophysicalModel), Cmu_ ( dimensioned::lookupOrAddToDict ( "Cmu", coeffDict_, 0.09 ) ), C1_ ( dimensioned::lookupOrAddToDict ( "C1", coeffDict_, 1.44 ) ), C2_ ( dimensioned::lookupOrAddToDict ( "C2", coeffDict_, 1.92 ) ), sigmak_ ( dimensioned::lookupOrAddToDict ( "sigmak", coeffDict_, 1.0 ) ), sigmaEps_ ( dimensioned::lookupOrAddToDict ( "sigmaEps", coeffDict_, 1.3 ) ), Prt_ ( dimensioned::lookupOrAddToDict ( "Prt", coeffDict_, 1.0 ) ), 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_ ), mut_ ( IOobject ( "mut", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), Cmu_*rho_*sqr(k_)/(epsilon_ + epsilonSmall_) ), alphat_ ( IOobject ( "alphat", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::AUTO_WRITE ), autoCreateAlphat("alphat", mesh_) ) { mut_ = Cmu_*rho_*sqr(k_)/(epsilon_ + epsilonSmall_); mut_.correctBoundaryConditions(); alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); printCoeffs(); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // tmp PDRkEpsilon::R() const { return tmp ( new volSymmTensorField ( IOobject ( "R", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), ((2.0/3.0)*I)*k_ - (mut_/rho_)*dev(twoSymm(fvc::grad(U_))), k_.boundaryField().types() ) ); } tmp PDRkEpsilon::devRhoReff() const { return tmp ( new volSymmTensorField ( IOobject ( "devRhoReff", runTime_.timeName(), mesh_, IOobject::NO_READ, IOobject::NO_WRITE ), -muEff()*dev(twoSymm(fvc::grad(U_))) ) ); } tmp PDRkEpsilon::divDevRhoReff(volVectorField& U) const { return ( - fvm::laplacian(muEff(), U) - fvc::div(muEff()*dev2(T(fvc::grad(U)))) ); } bool PDRkEpsilon::read() { if (RASModel::read()) { Cmu_.readIfPresent(coeffDict_); C1_.readIfPresent(coeffDict_); C2_.readIfPresent(coeffDict_); sigmak_.readIfPresent(coeffDict()); sigmaEps_.readIfPresent(coeffDict()); Prt_.readIfPresent(coeffDict()); return true; } else { return false; } } void PDRkEpsilon::correct() { if (!turbulence_) { // Re-calculate viscosity mut_ = rho_*Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_); mut_.correctBoundaryConditions(); // Re-calculate thermal diffusivity alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); return; } RASModel::correct(); volScalarField divU = fvc::div(phi_/fvc::interpolate(rho_)); if (mesh_.moving()) { divU += fvc::div(mesh_.phi()); } tmp tgradU = fvc::grad(U_); volScalarField G = 2*mut_*(tgradU() && dev(symm(tgradU()))); tgradU.clear(); // Update espsilon and G at the wall epsilon_.boundaryField().updateCoeffs(); // Add the blockage generation term so that it is included consistently // in both the k and epsilon equations const volScalarField& betav = U_.db().lookupObject("betav"); const PDRDragModel& drag = U_.db().lookupObject("PDRDragModel"); volScalarField GR = drag.Gk(); // Dissipation equation tmp epsEqn ( betav*fvm::ddt(rho_, epsilon_) + fvm::div(phi_, epsilon_) - fvm::laplacian(DepsilonEff(), epsilon_) == C1_*(betav*G + GR)*epsilon_/k_ - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_) - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_) ); epsEqn().relax(); // No longer needed: matrix completes at the point of solution // HJ, 17/Apr/2012 // epsEqn().completeAssembly(); solve(epsEqn); bound(epsilon_, epsilon0_); // Turbulent kinetic energy equation tmp kEqn ( betav*fvm::ddt(rho_, k_) + fvm::div(phi_, k_) - fvm::laplacian(DkEff(), k_) == betav*G + GR - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_) - fvm::Sp(betav*rho_*epsilon_/k_, k_) ); kEqn().relax(); solve(kEqn); bound(k_, k0_); // Re-calculate viscosity mut_ = rho_*Cmu_*sqr(k_)/epsilon_; mut_.correctBoundaryConditions(); // Re-calculate thermal diffusivity alphat_ = mut_/Prt_; alphat_.correctBoundaryConditions(); } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // } // End namespace RASModels } // End namespace compressible } // End namespace Foam // ************************************************************************* //