Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / LES / SpalartAllmarasIDDES / IDDESDelta / IDDESDelta.C
blobf1e008eab1ff54db60503f2399274300ffa2e65f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-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 "IDDESDelta.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "wallDistReflection.H"
29 #include "wallDist.H"
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace Foam
35     defineTypeNameAndDebug(IDDESDelta, 0);
36     addToRunTimeSelectionTable(LESdelta, IDDESDelta, dictionary);
40 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
42 void Foam::IDDESDelta::calcDelta()
44     label nD = mesh().nGeometricD();
46     const volScalarField& hmax = hmax_();
48     // initialise wallNorm
49     wallDistReflection wallNorm(mesh());
51     const volVectorField& n = wallNorm.n();
53     tmp<volScalarField> faceToFacenMax
54     (
55         new volScalarField
56         (
57             IOobject
58             (
59                 "faceToFaceMax",
60                 mesh().time().timeName(),
61                 mesh(),
62                 IOobject::NO_READ,
63                 IOobject::NO_WRITE
64             ),
65             mesh(),
66             dimensionedScalar("zrero", dimLength, 0.0)
67         )
68     );
70     const cellList& cells = mesh().cells();
72     forAll(cells,cellI)
73     {
74         scalar deltaMaxTmp = 0.0;
75         const labelList& cFaces = mesh().cells()[cellI];
76         const point& faceCentre = mesh().faceCentres()[cFaces[0]];
77         const vector nCell = n[cellI];
78         forAll(cFaces, cFaceI)
79         {
80             label faceI = cFaces[cFaceI];
81             const point& faceCentreTwo = mesh().faceCentres()[faceI];
82             scalar tmp = (faceCentre - faceCentreTwo) & nCell;
83             if (tmp > deltaMaxTmp)
84             {
85                 deltaMaxTmp = tmp;
86             }
87         }
88         faceToFacenMax()[cellI] = deltaMaxTmp;
89     }
91     if (nD == 3)
92     {
93         delta_.internalField() =
94             deltaCoeff_
95            *min
96             (
97                 max
98                 (
99                     max(cw_*wallDist(mesh()).y(), cw_*hmax),
100                     faceToFacenMax()
101                 ),
102                 hmax
103             );
104     }
105     else if (nD == 2)
106     {
107         WarningIn("IDDESDelta::calcDelta()")
108             << "Case is 2D, LES is not strictly applicable\n"
109             << endl;
111         delta_.internalField() =
112             deltaCoeff_
113            *min
114             (
115                 max(max(cw_*wallDist(mesh()).y(), cw_*hmax), faceToFacenMax()),
116                 hmax
117             );
118     }
119     else
120     {
121         FatalErrorIn("IDDESDelta::calcDelta()")
122             << "Case is not 3D or 2D, LES is not strictly applicable"
123             << exit(FatalError);
124     }
128 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
130 Foam::IDDESDelta::IDDESDelta
132     const word& name,
133     const fvMesh& mesh,
134     const dictionary& dd
137     LESdelta(name, mesh),
138     hmax_(LESdelta::New("hmax", mesh, dd.parent())),
139     deltaCoeff_
140     (
141         readScalar(dd.subDict(type()+"Coeffs").lookup("deltaCoeff"))
142     ),
143     cw_(0.15)
145     dd.subDict(type() + "Coeffs").readIfPresent("cw", cw_);
146     calcDelta();
150 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
152 void Foam::IDDESDelta::read(const dictionary& dd)
154     dd.subDict(type() + "Coeffs").lookup("deltaCoeff") >> deltaCoeff_;
155     calcDelta();
159 void Foam::IDDESDelta::correct()
161     if (mesh_.changing())
162     {
163         calcDelta();
164         hmax_().correct();
165     }
169 // ************************************************************************* //