initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / LES / vanDriestDelta / vanDriestDelta.C
blobe366ac159fd1fd0bcd357d4501ca5b8c1d3a7fb7
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 "vanDriestDelta.H"
28 #include "LESModel.H"
29 #include "wallFvPatch.H"
30 #include "wallDistData.H"
31 #include "wallPointYPlus.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace Foam
38 namespace compressible
40 namespace LESModels
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 defineTypeNameAndDebug(vanDriestDelta, 0);
46 addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
48 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
50 void vanDriestDelta::calcDelta()
52     const LESModel& lesModel = mesh_.lookupObject<LESModel>("LESProperties");
54     const volVectorField& U = lesModel.U();
55     const volScalarField& rho = lesModel.rho();
56     const volScalarField& mu = lesModel.mu();
57     tmp<volScalarField> muSgs = lesModel.muSgs();
59     volScalarField ystar
60     (
61         IOobject
62         (
63             "ystar",
64             mesh_.time().constant(),
65             mesh_
66         ),
67         mesh_,
68         dimensionedScalar("ystar", dimLength, GREAT)
69     );
71     const fvPatchList& patches = mesh_.boundary();
72     forAll(patches, patchi)
73     {
74         if (isType<wallFvPatch>(patches[patchi]))
75         {
76             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
77             const scalarField& rhow = rho.boundaryField()[patchi];
78             const scalarField& muw = mu.boundaryField()[patchi];
79             const scalarField& muSgsw = muSgs().boundaryField()[patchi];
81             ystar.boundaryField()[patchi] =
82                 muw/(rhow*sqrt((muw + muSgsw)*mag(Uw.snGrad())/rhow + VSMALL));
83         }
84     }
86     wallPointYPlus::yPlusCutOff = 500;
87     wallDistData<wallPointYPlus> y(mesh_, ystar);
89     delta_ = min
90     (
91         static_cast<const volScalarField&>(geometricDelta_()),
92         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
93     );
97 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
99 vanDriestDelta::vanDriestDelta
101     const word& name,
102     const fvMesh& mesh,
103     const dictionary& dd
106     LESdelta(name, mesh),
107     geometricDelta_
108     (
109         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
110     ),
111     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
112     Aplus_
113     (
114         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
115     ),
116     Cdelta_
117     (
118         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
119     ),
120     calcInterval_
121     (
122         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
123     )
125     delta_ = geometricDelta_();
129 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
131 void vanDriestDelta::read(const dictionary& d)
133     const dictionary& dd(d.subDict(type() + "Coeffs"));
135     geometricDelta_().read(dd);
136     d.readIfPresent<scalar>("kappa", kappa_);
137     dd.readIfPresent<scalar>("Aplus", Aplus_);
138     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
139     dd.readIfPresent<label>("calcInterval", calcInterval_);
140     calcDelta();
144 void vanDriestDelta::correct()
146     if (mesh().time().timeIndex() % calcInterval_ == 0)
147     {
148         geometricDelta_().correct();
149         calcDelta();
150     }
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 } // End namespace LESModels
157 } // End namespace compressible
158 } // End namespace Foam
160 // ************************************************************************* //