1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "DeardorffDiffStress.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace incompressible
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(DeardorffDiffStress, 0);
42 addToRunTimeSelectionTable(LESModel, DeardorffDiffStress, dictionary);
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 void DeardorffDiffStress::updateSubGridScaleFields(const volScalarField& K)
48 nuSgs_ = ck_*sqrt(K)*delta();
49 nuSgs_.correctBoundaryConditions();
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 DeardorffDiffStress::DeardorffDiffStress
57 const volVectorField& U,
58 const surfaceScalarField& phi,
59 transportModel& transport
62 LESModel(typeName, U, phi, transport),
63 GenSGSStress(U, phi, transport),
67 dimensioned<scalar>::lookupOrAddToDict
76 dimensioned<scalar>::lookupOrAddToDict
84 updateSubGridScaleFields(0.5*tr(B_));
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 void DeardorffDiffStress::correct(const tmp<volTensorField>& tgradU)
94 const volTensorField& gradU = tgradU();
96 GenSGSStress::correct(gradU);
98 volSymmTensorField D = symm(gradU);
100 volSymmTensorField P = -twoSymm(B_ & gradU);
102 volScalarField K = 0.5*tr(B_);
103 volScalarField Epsilon = 2*nuEff()*magSqr(D);
105 fvSymmTensorMatrix BEqn
108 + fvm::div(phi(), B_)
109 - fvm::laplacian(DBEff(), B_)
110 + fvm::Sp(cm_*sqrt(K)/delta(), B_)
114 - (2*ce_ - 0.667*cm_)*I*Epsilon
120 // Bounding the component kinetic energies
124 B_[celli].component(symmTensor::XX) =
125 max(B_[celli].component(symmTensor::XX), k0().value());
126 B_[celli].component(symmTensor::YY) =
127 max(B_[celli].component(symmTensor::YY), k0().value());
128 B_[celli].component(symmTensor::ZZ) =
129 max(B_[celli].component(symmTensor::ZZ), k0().value());
135 updateSubGridScaleFields(K);
139 bool DeardorffDiffStress::read()
141 if (GenSGSStress::read())
143 ck_.readIfPresent(coeffDict());
144 cm_.readIfPresent(coeffDict());
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 } // End namespace LESModels
158 } // End namespace incompressible
159 } // End namespace Foam
161 // ************************************************************************* //