initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / RAS / backwardsCompatibility / wallFunctions / backwardsCompatibilityWallFunctions.C
blob7096314fe34aaddbeea31b6e88cab49a70719bf6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-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 "backwardsCompatibilityWallFunctions.H"
29 #include "calculatedFvPatchField.H"
30 #include "alphatWallFunctionFvPatchScalarField.H"
31 #include "mutWallFunctionFvPatchScalarField.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 (isType<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 (isType<wallFvPatch>(bm[patchI]))
147             {
148                 mutBoundaryTypes[patchI] =
149                     RASModels::mutWallFunctionFvPatchScalarField::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> autoCreateEpsilon
187     const word& fieldName,
188     const fvMesh& mesh
191     return
192         autoCreateWallFunctionField
193         <
194             scalar,
195             RASModels::epsilonWallFunctionFvPatchScalarField
196         >
197         (
198             fieldName,
199             mesh
200         );
204 tmp<volScalarField> autoCreateOmega
206     const word& fieldName,
207     const fvMesh& mesh
210     return
211         autoCreateWallFunctionField
212         <
213             scalar,
214             RASModels::omegaWallFunctionFvPatchScalarField
215         >
216         (
217             fieldName,
218             mesh
219         );
223 tmp<volScalarField> autoCreateK
225     const word& fieldName,
226     const fvMesh& mesh
229     return
230         autoCreateWallFunctionField
231         <
232             scalar,
233             RASModels::kqRWallFunctionFvPatchField<scalar>
234         >
235         (
236             fieldName,
237             mesh
238         );
242 tmp<volScalarField> autoCreateQ
244     const word& fieldName,
245     const fvMesh& mesh
248     return
249         autoCreateWallFunctionField
250         <
251             scalar,
252             RASModels::kqRWallFunctionFvPatchField<scalar>
253         >
254         (
255             fieldName,
256             mesh
257         );
261 tmp<volSymmTensorField> autoCreateR
263     const word& fieldName,
264     const fvMesh& mesh
267     return
268         autoCreateWallFunctionField
269         <
270             symmTensor,
271             RASModels::kqRWallFunctionFvPatchField<symmTensor>
272         >
273         (
274             fieldName,
275             mesh
276         );
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 } // End namespace compressible
283 } // End namespace Foam
285 // ************************************************************************* //