Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / compressible / RAS / backwardsCompatibility / wallFunctions / backwardsCompatibilityWallFunctions.C
blobbdf8d01de3ff3813f17c0a85f4f0b9fb97b88a59
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 "backwardsCompatibilityWallFunctions.H"
28 #include "calculatedFvPatchField.H"
29 #include "alphatWallFunctionFvPatchScalarField.H"
30 #include "mutkWallFunctionFvPatchScalarField.H"
31 #include "mutLowReWallFunctionFvPatchScalarField.H"
32 #include "epsilonWallFunctionFvPatchScalarField.H"
33 #include "kqRWallFunctionFvPatchField.H"
34 #include "omegaWallFunctionFvPatchScalarField.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace Foam
40 namespace compressible
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 tmp<volScalarField> autoCreateAlphat
47     const word& fieldName,
48     const fvMesh& mesh
51     IOobject alphatHeader
52     (
53         fieldName,
54         mesh.time().timeName(),
55         mesh,
56         IOobject::MUST_READ,
57         IOobject::NO_WRITE,
58         false
59     );
61     if (alphatHeader.headerOk())
62     {
63         return tmp<volScalarField>(new volScalarField(alphatHeader, mesh));
64     }
65     else
66     {
67         Info<< "--> Creating " << fieldName
68             << " to employ run-time selectable wall functions" << endl;
70         const fvBoundaryMesh& bm = mesh.boundary();
72         wordList alphatBoundaryTypes(bm.size());
74         forAll(bm, patchI)
75         {
76             if (isA<wallFvPatch>(bm[patchI]))
77             {
78                 alphatBoundaryTypes[patchI] =
79                     RASModels::alphatWallFunctionFvPatchScalarField::typeName;
80             }
81             else
82             {
83                 alphatBoundaryTypes[patchI] =
84                     calculatedFvPatchField<scalar>::typeName;
85             }
86         }
88         tmp<volScalarField> alphat
89         (
90             new volScalarField
91             (
92                 IOobject
93                 (
94                     fieldName,
95                     mesh.time().timeName(),
96                     mesh,
97                     IOobject::NO_READ,
98                     IOobject::NO_WRITE,
99                     false
100                 ),
101                 mesh,
102                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
103                 alphatBoundaryTypes
104             )
105         );
107         Info<< "    Writing new " << fieldName << endl;
108         alphat().write();
110         return alphat;
111     }
115 tmp<volScalarField> autoCreateMut
117     const word& fieldName,
118     const fvMesh& mesh
121     IOobject mutHeader
122     (
123         fieldName,
124         mesh.time().timeName(),
125         mesh,
126         IOobject::MUST_READ,
127         IOobject::NO_WRITE,
128         false
129     );
131     if (mutHeader.headerOk())
132     {
133         return tmp<volScalarField>(new volScalarField(mutHeader, mesh));
134     }
135     else
136     {
137         Info<< "--> Creating " << fieldName
138             << " to employ run-time selectable wall functions" << endl;
140         const fvBoundaryMesh& bm = mesh.boundary();
142         wordList mutBoundaryTypes(bm.size());
144         forAll(bm, patchI)
145         {
146             if (isA<wallFvPatch>(bm[patchI]))
147             {
148                 mutBoundaryTypes[patchI] =
149                     RASModels::mutkWallFunctionFvPatchScalarField::typeName;
150             }
151             else
152             {
153                 mutBoundaryTypes[patchI] =
154                     calculatedFvPatchField<scalar>::typeName;
155             }
156         }
158         tmp<volScalarField> mut
159         (
160             new volScalarField
161             (
162                 IOobject
163                 (
164                     fieldName,
165                     mesh.time().timeName(),
166                     mesh,
167                     IOobject::NO_READ,
168                     IOobject::NO_WRITE,
169                     false
170                 ),
171                 mesh,
172                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
173                 mutBoundaryTypes
174             )
175         );
177         Info<< "    Writing new " << fieldName << endl;
178         mut().write();
180         return mut;
181     }
185 tmp<volScalarField> autoCreateLowReMut
187     const word& fieldName,
188     const fvMesh& mesh
191     IOobject mutHeader
192     (
193         fieldName,
194         mesh.time().timeName(),
195         mesh,
196         IOobject::MUST_READ,
197         IOobject::NO_WRITE,
198         false
199     );
201     if (mutHeader.headerOk())
202     {
203         return tmp<volScalarField>(new volScalarField(mutHeader, mesh));
204     }
205     else
206     {
207         Info<< "--> Creating " << fieldName
208             << " to employ run-time selectable wall functions" << endl;
210         const fvBoundaryMesh& bm = mesh.boundary();
212         wordList mutBoundaryTypes(bm.size());
214         forAll(bm, patchI)
215         {
216             if (isA<wallFvPatch>(bm[patchI]))
217             {
218                 mutBoundaryTypes[patchI] =
219                     RASModels::mutLowReWallFunctionFvPatchScalarField::typeName;
220             }
221             else
222             {
223                 mutBoundaryTypes[patchI] =
224                     calculatedFvPatchField<scalar>::typeName;
225             }
226         }
228         tmp<volScalarField> mut
229         (
230             new volScalarField
231             (
232                 IOobject
233                 (
234                     fieldName,
235                     mesh.time().timeName(),
236                     mesh,
237                     IOobject::NO_READ,
238                     IOobject::NO_WRITE,
239                     false
240                 ),
241                 mesh,
242                 dimensionedScalar("zero", dimDensity*dimArea/dimTime, 0.0),
243                 mutBoundaryTypes
244             )
245         );
247         Info<< "    Writing new " << fieldName << endl;
248         mut().write();
250         return mut;
251     }
255 tmp<volScalarField> autoCreateEpsilon
257     const word& fieldName,
258     const fvMesh& mesh
261     return
262         autoCreateWallFunctionField
263         <
264             scalar,
265             RASModels::epsilonWallFunctionFvPatchScalarField
266         >
267         (
268             fieldName,
269             mesh
270         );
274 tmp<volScalarField> autoCreateOmega
276     const word& fieldName,
277     const fvMesh& mesh
280     return
281         autoCreateWallFunctionField
282         <
283             scalar,
284             RASModels::omegaWallFunctionFvPatchScalarField
285         >
286         (
287             fieldName,
288             mesh
289         );
293 tmp<volScalarField> autoCreateK
295     const word& fieldName,
296     const fvMesh& mesh
299     return
300         autoCreateWallFunctionField
301         <
302             scalar,
303             RASModels::kqRWallFunctionFvPatchField<scalar>
304         >
305         (
306             fieldName,
307             mesh
308         );
312 tmp<volScalarField> autoCreateQ
314     const word& fieldName,
315     const fvMesh& mesh
318     return
319         autoCreateWallFunctionField
320         <
321             scalar,
322             RASModels::kqRWallFunctionFvPatchField<scalar>
323         >
324         (
325             fieldName,
326             mesh
327         );
331 tmp<volSymmTensorField> autoCreateR
333     const word& fieldName,
334     const fvMesh& mesh
337     return
338         autoCreateWallFunctionField
339         <
340             symmTensor,
341             RASModels::kqRWallFunctionFvPatchField<symmTensor>
342         >
343         (
344             fieldName,
345             mesh
346         );
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 } // End namespace compressible
353 } // End namespace Foam
355 // ************************************************************************* //