ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / incompressible / LES / vanDriestDelta / vanDriestDelta.C
blobf6c77ad9838fd71929b691d2f3beb4fed7c7b41b
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
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
19     for more details.
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"
27 #include "LESModel.H"
28 #include "wallFvPatch.H"
29 #include "wallDistData.H"
30 #include "wallPointYPlus.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
37 namespace incompressible
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& nu = lesModel.nu();
55     tmp<volScalarField> nuSgs = lesModel.nuSgs();
57     volScalarField ystar
58     (
59         IOobject
60         (
61             "ystar",
62             mesh_.time().constant(),
63             mesh_
64         ),
65         mesh_,
66         dimensionedScalar("ystar", dimLength, GREAT)
67     );
69     const fvPatchList& patches = mesh_.boundary();
70     forAll(patches, patchi)
71     {
72         if (isA<wallFvPatch>(patches[patchi]))
73         {
74             const fvPatchVectorField& Uw = U.boundaryField()[patchi];
75             const scalarField& nuw = nu.boundaryField()[patchi];
76             const scalarField& nuSgsw = nuSgs().boundaryField()[patchi];
78             ystar.boundaryField()[patchi] =
79                 nuw/sqrt((nuw + nuSgsw)*mag(Uw.snGrad()) + VSMALL);
80         }
81     }
83     wallPointYPlus::yPlusCutOff = 500;
84     wallDistData<wallPointYPlus> y(mesh_, ystar);
86     delta_ = min
87     (
88         static_cast<const volScalarField&>(geometricDelta_()),
89         (kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
90     );
94 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
96 vanDriestDelta::vanDriestDelta
98     const word& name,
99     const fvMesh& mesh,
100     const dictionary& dd
103     LESdelta(name, mesh),
104     geometricDelta_
105     (
106         LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
107     ),
108     kappa_(dd.lookupOrDefault<scalar>("kappa", 0.41)),
109     Aplus_
110     (
111         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
112     ),
113     Cdelta_
114     (
115         dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
116     ),
117     calcInterval_
118     (
119         dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
120     )
122     delta_ = geometricDelta_();
126 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
128 void vanDriestDelta::read(const dictionary& d)
130     const dictionary& dd(d.subDict(type() + "Coeffs"));
132     geometricDelta_().read(dd);
133     d.readIfPresent<scalar>("kappa", kappa_);
134     dd.readIfPresent<scalar>("Aplus", Aplus_);
135     dd.readIfPresent<scalar>("Cdelta", Cdelta_);
136     dd.readIfPresent<label>("calcInterval", calcInterval_);
137     calcDelta();
141 void vanDriestDelta::correct()
143     if (mesh().time().timeIndex() % calcInterval_ == 0)
144     {
145         geometricDelta_().correct();
146         calcDelta();
147     }
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 } // End namespace LESModels
154 } // End namespace incompressible
155 } // End namespace Foam
157 // ************************************************************************* //