/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2013-2017 OpenFOAM Foundation Copyright (C) 2019-2021 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM 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. OpenFOAM 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 OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "LESModel.H" // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * // template void Foam::LESModel::printCoeffs(const word& type) { if (printCoeffs_) { Info<< coeffDict_.dictName() << coeffDict_ << endl; } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // template Foam::LESModel::LESModel ( const word& type, const alphaField& alpha, const rhoField& rho, const volVectorField& U, const surfaceScalarField& alphaRhoPhi, const surfaceScalarField& phi, const transportModel& transport, const word& propertiesName ) : BasicTurbulenceModel ( type, alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName ), LESDict_(this->subOrEmptyDict("LES")), turbulence_(LESDict_.getOrDefault("turbulence", true)), printCoeffs_(LESDict_.getOrDefault("printCoeffs", false)), coeffDict_(LESDict_.optionalSubDict(type + "Coeffs")), Ce_ ( dimensioned::getOrAddToDict ( "Ce", LESDict_, 1.048 ) ), kMin_ ( dimensioned::getOrAddToDict ( "kMin", LESDict_, sqr(dimVelocity), SMALL ) ), epsilonMin_ ( dimensioned::getOrAddToDict ( "epsilonMin", LESDict_, kMin_.dimensions()/dimTime, SMALL ) ), omegaMin_ ( dimensioned::getOrAddToDict ( "omegaMin", LESDict_, dimless/dimTime, SMALL ) ), delta_ ( LESdelta::New ( IOobject::groupName("delta", alphaRhoPhi.group()), *this, LESDict_ ) ) { // Force the construction of the mesh deltaCoeffs which may be needed // for the construction of the derived models and BCs this->mesh_.deltaCoeffs(); } // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * // template Foam::autoPtr> Foam::LESModel::New ( const alphaField& alpha, const rhoField& rho, const volVectorField& U, const surfaceScalarField& alphaRhoPhi, const surfaceScalarField& phi, const transportModel& transport, const word& propertiesName ) { const IOdictionary modelDict ( IOobject ( IOobject::groupName(propertiesName, alphaRhoPhi.group()), U.time().constant(), U.db(), IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER ) ); const dictionary& dict = modelDict.subDict("LES"); const word modelType ( // "LESModel" for v2006 and earlier dict.getCompat("model", {{"LESModel", -2006}}) ); Info<< "Selecting LES turbulence model " << modelType << endl; auto* ctorPtr = dictionaryConstructorTable(modelType); if (!ctorPtr) { FatalIOErrorInLookup ( dict, "LES model", modelType, *dictionaryConstructorTablePtr_ ) << exit(FatalIOError); } return autoPtr ( ctorPtr(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName) ); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // template bool Foam::LESModel::read() { if (BasicTurbulenceModel::read()) { LESDict_ <<= this->subDict("LES"); LESDict_.readEntry("turbulence", turbulence_); coeffDict_ <<= LESDict_.optionalSubDict(type() + "Coeffs"); delta_().read(LESDict_); Ce_.readIfPresent(LESDict_); kMin_.readIfPresent(LESDict_); return true; } return false; } template Foam::tmp Foam::LESModel::epsilon() const { return tmp::New ( IOobject ( IOobject::groupName("epsilon", this->alphaRhoPhi_.group()), this->mesh_.time().timeName(), this->mesh_ ), this->Ce()*pow(this->k(), 1.5)/this->delta() ); } template Foam::tmp Foam::LESModel::omega() const { const scalar betaStar = 0.09; const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL); return tmp::New ( IOobject ( IOobject::groupName("omega", this->alphaRhoPhi_.group()), this->mesh_.time().timeName(), this->mesh_ ), this->epsilon()/(betaStar*(this->k() + k0)) ); } template void Foam::LESModel::correct() { delta_().correct(); BasicTurbulenceModel::correct(); } // ************************************************************************* //