Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / turbulenceModels / incompressible / LES / vanDriestDelta / vanDriestDelta.C
blobdd453bd9847850a51e79d9ffac7ad3a32f6c52eb
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "vanDriestDelta.H"
27 #include "LESModel.H"
28 #include "wallDistData.H"
29 #include "wallPointYPlus.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
36 namespace incompressible
38 namespace LESModels
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 defineTypeNameAndDebug(vanDriestDelta, 0);
44 addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
46 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
48 void vanDriestDelta::calcDelta()
50     const LESModel& lesModel = mesh_.lookupObject<LESModel>("LESProperties");
52     const volVectorField& U = lesModel.U();
53     const volScalarField& nu = lesModel.nu();
54     tmp<volScalarField> nuSgs = lesModel.nuSgs();
56     volScalarField ystar
57     (
58         IOobject
59         (
60             "ystar",
61             mesh_.time().constant(),
62             U.db()
63         ),
64         mesh_,
65         dimensionedScalar("ystar", dimLength, GREAT)
66     );
68     const fvPatchList& patches = mesh_.boundary();
69     forAll(patches, patchi)
70     {
71         if (patches[patchi].isWall())
72         {
73             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
74             const scalarField& nuw = nu.boundaryField()[patchi];
75             const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
77             ystar.boundaryField()[patchi] =
78                 nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
79         }
80     }
82     wallPointYPlus::yPlusCutOff = 500;
83     wallDistData<wallPointYPlus> y(mesh_, ystar);
85     delta_ = min
86     (
87         static_cast<const volScalarField&>(geometricDelta_()),
88         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
89     );
93 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
95 vanDriestDelta::vanDriestDelta
97     const word& name,
98     const fvMesh& mesh,
99     const dictionary& dd
102     LESdelta(name, mesh),
103     geometricDelta_
104     (
105         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
106     ),
107     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
108     Aplus_
109     (
110         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
111     ),
112     Cdelta_
113     (
114         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
115     ),
116     calcInterval_
117     (
118         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
119     )
121     delta_ = geometricDelta_();
125 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
127 void vanDriestDelta::read(const dictionary& d)
129     const dictionary& dd(d.subDict(type() + "Coeffs"));
131     geometricDelta_().read(dd);
132     d.readIfPresent<scalar>("kappa", kappa_);
133     dd.readIfPresent<scalar>("Aplus", Aplus_);
134     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
135     dd.readIfPresent<label>("calcInterval", calcInterval_);
136     calcDelta();
140 void vanDriestDelta::correct()
142     if (mesh().time().timeIndex() % calcInterval_ == 0)
143     {
144         geometricDelta_().correct();
145         calcDelta();
146     }
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 } // End namespace LESModels
153 } // End namespace incompressible
154 } // End namespace Foam
156 // ************************************************************************* //