1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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 "vanDriestDelta.H"
28 #include "wallFvPatch.H"
29 #include "wallDistData.H"
30 #include "wallPointYPlus.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace incompressible
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 defineTypeNameAndDebug(vanDriestDelta, 0);
45 addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 void vanDriestDelta::calcDelta()
51 const LESModel& lesModel = mesh_.lookupObject<LESModel>("LESProperties");
53 const volVectorField& U = lesModel.U();
54 const tmp<volScalarField> tnu = lesModel.nu();
55 const volScalarField& nu = tnu();
56 tmp<volScalarField> nuSgs = lesModel.nuSgs();
63 mesh_.time().constant(),
67 dimensionedScalar("ystar", dimLength, GREAT)
70 const fvPatchList& patches = mesh_.boundary();
71 forAll(patches, patchi)
73 if (isA<wallFvPatch>(patches[patchi]))
75 const fvPatchVectorField& Uw = U.boundaryField()[patchi];
76 const scalarField& nuw = nu.boundaryField()[patchi];
77 const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
79 ystar.boundaryField()[patchi] =
80 nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
84 wallPointYPlus::yPlusCutOff = 500;
85 wallDistData<wallPointYPlus> y(mesh_, ystar);
89 static_cast<const volScalarField&>(geometricDelta_()),
90 (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 vanDriestDelta::vanDriestDelta
104 LESdelta(name, mesh),
107 LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
109 kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
112 dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
116 dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
120 dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
123 delta_ = geometricDelta_();
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 void vanDriestDelta::read(const dictionary& d)
131 const dictionary& dd(d.subDict(type() + "Coeffs"));
133 geometricDelta_().read(dd);
134 d.readIfPresent<scalar>("kappa", kappa_);
135 dd.readIfPresent<scalar>("Aplus", Aplus_);
136 dd.readIfPresent<scalar>("Cdelta", Cdelta_);
137 dd.readIfPresent<label>("calcInterval", calcInterval_);
142 void vanDriestDelta::correct()
144 if (mesh().time().timeIndex() % calcInterval_ == 0)
146 geometricDelta_().correct();
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 } // End namespace LESModels
155 } // End namespace incompressible
156 } // End namespace Foam
158 // ************************************************************************* //