Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / backwardsCompatibility / wallFunctions / backwardsCompatibilityWallFunctions.C
blobdc2394057afa2f5c1f96c6f69b5cdcaddd24c1bc
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 "nutkWallFunctionFvPatchScalarField.H"
30 #include "nutLowReWallFunctionFvPatchScalarField.H"
31 #include "epsilonWallFunctionFvPatchScalarField.H"
32 #include "kqRWallFunctionFvPatchField.H"
33 #include "omegaWallFunctionFvPatchScalarField.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
39 namespace incompressible
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 tmp<volScalarField> autoCreateNut
46     const word& fieldName,
47     const fvMesh& mesh
50     IOobject nutHeader
51     (
52         fieldName,
53         mesh.time().timeName(),
54         mesh,
55         IOobject::MUST_READ,
56         IOobject::NO_WRITE,
57         false
58     );
60     if (nutHeader.headerOk())
61     {
62         return tmp<volScalarField>(new volScalarField(nutHeader, mesh));
63     }
64     else
65     {
66         Info<< "--> Creating " << fieldName
67             << " to employ run-time selectable wall functions" << endl;
69         const fvBoundaryMesh& bm = mesh.boundary();
71         wordList nutBoundaryTypes(bm.size());
73         forAll(bm, patchI)
74         {
75             if (isA<wallFvPatch>(bm[patchI]))
76             {
77                 nutBoundaryTypes[patchI] =
78                     RASModels::nutkWallFunctionFvPatchScalarField::typeName;
79             }
80             else
81             {
82                 nutBoundaryTypes[patchI] =
83                     calculatedFvPatchField<scalar>::typeName;
84             }
85         }
87         tmp<volScalarField> nut
88         (
89             new volScalarField
90             (
91                 IOobject
92                 (
93                     fieldName,
94                     mesh.time().timeName(),
95                     mesh,
96                     IOobject::NO_READ,
97                     IOobject::NO_WRITE,
98                     false
99                 ),
100                 mesh,
101                 dimensionedScalar("zero", dimArea/dimTime, 0.0),
102                 nutBoundaryTypes
103             )
104         );
106         Info<< "    Writing new " << fieldName << endl;
107         nut().write();
109         return nut;
110     }
114 tmp<volScalarField> autoCreateLowReNut
116     const word& fieldName,
117     const fvMesh& mesh
120     IOobject nutHeader
121     (
122         fieldName,
123         mesh.time().timeName(),
124         mesh,
125         IOobject::MUST_READ,
126         IOobject::NO_WRITE,
127         false
128     );
130     if (nutHeader.headerOk())
131     {
132         return tmp<volScalarField>(new volScalarField(nutHeader, mesh));
133     }
134     else
135     {
136         Info<< "--> Creating " << fieldName
137             << " to employ run-time selectable wall functions" << endl;
139         const fvBoundaryMesh& bm = mesh.boundary();
141         wordList nutBoundaryTypes(bm.size());
143         forAll(bm, patchI)
144         {
145             if (isA<wallFvPatch>(bm[patchI]))
146             {
147                 nutBoundaryTypes[patchI] =
148                     RASModels::nutLowReWallFunctionFvPatchScalarField::typeName;
149             }
150             else
151             {
152                 nutBoundaryTypes[patchI] =
153                     calculatedFvPatchField<scalar>::typeName;
154             }
155         }
157         tmp<volScalarField> nut
158         (
159             new volScalarField
160             (
161                 IOobject
162                 (
163                     fieldName,
164                     mesh.time().timeName(),
165                     mesh,
166                     IOobject::NO_READ,
167                     IOobject::NO_WRITE,
168                     false
169                 ),
170                 mesh,
171                 dimensionedScalar("zero", dimArea/dimTime, 0.0),
172                 nutBoundaryTypes
173             )
174         );
176         Info<< "    Writing new " << fieldName << endl;
177         nut().write();
179         return nut;
180     }
184 tmp<volScalarField> autoCreateEpsilon
186     const word& fieldName,
187     const fvMesh& mesh
190     return
191         autoCreateWallFunctionField
192         <
193             scalar,
194             RASModels::epsilonWallFunctionFvPatchScalarField
195         >
196         (
197             fieldName,
198             mesh
199         );
203 tmp<volScalarField> autoCreateOmega
205     const word& fieldName,
206     const fvMesh& mesh
209     return
210         autoCreateWallFunctionField
211         <
212             scalar,
213             RASModels::omegaWallFunctionFvPatchScalarField
214         >
215         (
216             fieldName,
217             mesh
218         );
222 tmp<volScalarField> autoCreateK
224     const word& fieldName,
225     const fvMesh& mesh
228     return
229         autoCreateWallFunctionField
230         <
231             scalar,
232             RASModels::kqRWallFunctionFvPatchField<scalar>
233         >
234         (
235             fieldName,
236             mesh
237         );
241 tmp<volScalarField> autoCreateQ
243     const word& fieldName,
244     const fvMesh& mesh
247     return
248         autoCreateWallFunctionField
249         <
250             scalar,
251             RASModels::kqRWallFunctionFvPatchField<scalar>
252         >
253         (
254             fieldName,
255             mesh
256         );
260 tmp<volSymmTensorField> autoCreateR
262     const word& fieldName,
263     const fvMesh& mesh
266     return
267         autoCreateWallFunctionField
268         <
269             symmTensor,
270             RASModels::kqRWallFunctionFvPatchField<symmTensor>
271         >
272         (
273             fieldName,
274             mesh
275         );
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 } // End namespace incompressible
282 } // End namespace Foam
284 // ************************************************************************* //