/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2019-2023 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 "dynamicLagrangian.H"
#include "fvOptions.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace LESModels
{
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
template
void dynamicLagrangian::correctNut
(
const tmp& gradU
)
{
this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(devSymm(gradU));
this->nut_.correctBoundaryConditions();
fv::options::New(this->mesh_).correct(this->nut_);
BasicTurbulenceModel::correctNut();
}
template
void dynamicLagrangian::correctNut()
{
correctNut(fvc::grad(this->U_));
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
dynamicLagrangian::dynamicLagrangian
(
const alphaField& alpha,
const rhoField& rho,
const volVectorField& U,
const surfaceScalarField& alphaRhoPhi,
const surfaceScalarField& phi,
const transportModel& transport,
const word& propertiesName,
const word& type
)
:
LESeddyViscosity
(
type,
alpha,
rho,
U,
alphaRhoPhi,
phi,
transport,
propertiesName
),
flm_
(
IOobject
(
IOobject::groupName("flm", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
this->mesh_
),
fmm_
(
IOobject
(
IOobject::groupName("fmm", this->alphaRhoPhi_.group()),
this->runTime_.timeName(),
this->mesh_,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
this->mesh_
),
theta_
(
dimensioned::getOrAddToDict
(
"theta",
this->coeffDict_,
1.5
)
),
simpleFilter_(U.mesh()),
filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
filter_(filterPtr_()),
flm0_("flm0", flm_.dimensions(), Zero),
fmm0_("fmm0", fmm_.dimensions(), VSMALL)
{
if (type == typeName)
{
this->printCoeffs(type);
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
bool dynamicLagrangian::read()
{
if (LESeddyViscosity::read())
{
filter_.read(this->coeffDict());
theta_.readIfPresent(this->coeffDict());
return true;
}
return false;
}
template
void dynamicLagrangian::correct()
{
if (!this->turbulence_)
{
return;
}
// Local references
const alphaField& alpha = this->alpha_;
const rhoField& rho = this->rho_;
const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
const volVectorField& U = this->U_;
fv::options& fvOptions(fv::options::New(this->mesh_));
LESeddyViscosity::correct();
tmp tgradU(fvc::grad(U));
const volTensorField& gradU = tgradU();
volSymmTensorField S(devSymm(gradU));
volScalarField magS(mag(S));
volVectorField Uf(filter_(U));
volSymmTensorField Sf(devSymm(fvc::grad(Uf)));
volScalarField magSf(mag(Sf));
volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
volSymmTensorField M
(
2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
);
volScalarField invT
(
alpha*rho*(1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
);
volScalarField LM(L && M);
fvScalarMatrix flmEqn
(
fvm::ddt(alpha, rho, flm_)
+ fvm::div(alphaRhoPhi, flm_)
==
invT*LM
- fvm::Sp(invT, flm_)
+ fvOptions(alpha, rho, flm_)
);
flmEqn.relax();
fvOptions.constrain(flmEqn);
flmEqn.solve();
fvOptions.correct(flm_);
bound(flm_, flm0_);
volScalarField MM(M && M);
fvScalarMatrix fmmEqn
(
fvm::ddt(alpha, rho, fmm_)
+ fvm::div(alphaRhoPhi, fmm_)
==
invT*MM
- fvm::Sp(invT, fmm_)
+ fvOptions(alpha, rho, fmm_)
);
fmmEqn.relax();
fvOptions.constrain(fmmEqn);
fmmEqn.solve();
fvOptions.correct(fmm_);
bound(fmm_, fmm0_);
correctNut(gradU);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace LESModels
} // End namespace Foam
// ************************************************************************* //