changed isA<wallPolyPatch> & isA<wallFvPatch> almost everywhere
[openfoam-extend-OpenFOAM-1.6-ext.git] / src / turbulenceModels / compressible / LES / vanDriestDelta / vanDriestDelta.C
blob39f8a20d9a1e9b6fc64d8551bb1828da3a2aaf93
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "wallDistData.H"
30 #include "wallPointYPlus.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
37 namespace compressible
39 namespace LESModels
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 volScalarField& rho = lesModel.rho();
55     const volScalarField& mu = lesModel.mu();
56     tmp<volScalarField> muSgs = lesModel.muSgs();
58     volScalarField ystar
59     (
60         IOobject
61         (
62             "ystar",
63             mesh_.time().constant(),
64             mesh_
65         ),
66         mesh_,
67         dimensionedScalar("ystar", dimLength, GREAT)
68     );
70     const fvPatchList& patches = mesh_.boundary();
71     forAll(patches, patchi)
72     {
73         if (patches[patchi].isWall())
74         {
75             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
76             const scalarField& rhow = rho.boundaryField()[patchi];
77             const scalarField& muw = mu.boundaryField()[patchi];
78             const scalarField& muSgsw = muSgs().boundaryField()[patchi];
80             ystar.boundaryField()[patchi] =
81                 muw/(rhow*sqrt((muw + muSgsw)*mag(Uw.snGrad())/rhow + VSMALL));
82         }
83     }
85     wallPointYPlus::yPlusCutOff = 500;
86     wallDistData<wallPointYPlus> y(mesh_, ystar);
88     delta_ = min
89     (
90         static_cast<const volScalarField&>(geometricDelta_()),
91         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
92     );
96 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
98 vanDriestDelta::vanDriestDelta
100     const word& name,
101     const fvMesh& mesh,
102     const dictionary& dd
105     LESdelta(name, mesh),
106     geometricDelta_
107     (
108         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
109     ),
110     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
111     Aplus_
112     (
113         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
114     ),
115     Cdelta_
116     (
117         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
118     ),
119     calcInterval_
120     (
121         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
122     )
124     delta_ = geometricDelta_();
128 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
130 void vanDriestDelta::read(const dictionary& d)
132     const dictionary& dd(d.subDict(type() + "Coeffs"));
134     geometricDelta_().read(dd);
135     d.readIfPresent<scalar>("kappa", kappa_);
136     dd.readIfPresent<scalar>("Aplus", Aplus_);
137     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
138     dd.readIfPresent<label>("calcInterval", calcInterval_);
139     calcDelta();
143 void vanDriestDelta::correct()
145     if (mesh().time().timeIndex() % calcInterval_ == 0)
146     {
147         geometricDelta_().correct();
148         calcDelta();
149     }
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 } // End namespace LESModels
156 } // End namespace compressible
157 } // End namespace Foam
159 // ************************************************************************* //