1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "DeardorffDiffStress.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace incompressible
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(DeardorffDiffStress, 0);
41 addToRunTimeSelectionTable(LESModel, DeardorffDiffStress, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void DeardorffDiffStress::updateSubGridScaleFields(const volScalarField& K)
47 nuSgs_ = ck_*sqrt(K)*delta();
48 nuSgs_.correctBoundaryConditions();
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 DeardorffDiffStress::DeardorffDiffStress
56 const volVectorField& U,
57 const surfaceScalarField& phi,
58 transportModel& transport
61 LESModel(typeName, U, phi, transport),
62 GenSGSStress(U, phi, transport),
66 dimensioned<scalar>::lookupOrAddToDict
75 dimensioned<scalar>::lookupOrAddToDict
83 updateSubGridScaleFields(0.5*tr(B_));
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 void DeardorffDiffStress::correct(const tmp<volTensorField>& tgradU)
93 const volTensorField& gradU = tgradU();
95 GenSGSStress::correct(gradU);
97 volSymmTensorField D = symm(gradU);
99 volSymmTensorField P = -twoSymm(B_ & gradU);
101 volScalarField K = 0.5*tr(B_);
102 volScalarField Epsilon = 2*nuEff()*magSqr(D);
104 fvSymmTensorMatrix BEqn
107 + fvm::div(phi(), B_)
108 - fvm::laplacian(DBEff(), B_)
109 + fvm::Sp(cm_*sqrt(K)/delta(), B_)
113 - (2*ce_ - 0.667*cm_)*I*Epsilon
119 // Bounding the component kinetic energies
123 B_[celli].component(symmTensor::XX) =
124 max(B_[celli].component(symmTensor::XX), k0().value());
125 B_[celli].component(symmTensor::YY) =
126 max(B_[celli].component(symmTensor::YY), k0().value());
127 B_[celli].component(symmTensor::ZZ) =
128 max(B_[celli].component(symmTensor::ZZ), k0().value());
134 updateSubGridScaleFields(K);
138 bool DeardorffDiffStress::read()
140 if (GenSGSStress::read())
142 ck_.readIfPresent(coeffDict());
143 cm_.readIfPresent(coeffDict());
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 } // End namespace LESModels
157 } // End namespace incompressible
158 } // End namespace Foam
160 // ************************************************************************* //