Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / postProcessing / miscellaneous / execFlowFunctionObjects / execFlowFunctionObjects.C
blob2ead0ee53f401db2207ab2a0b4a366f72be155d4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Global
25     execFlowFunctionObjects
27 Description
28     Execute the set of functionObjects specified in the selected dictionary
29     (which defaults to system/controlDict) for the selected set of times.
31     The flow (p-U) and optionally turbulence fields are available for the
32     function objects to operate on allowing forces and other related properties
33     to be calculated in addition to cutting planes etc.
35 \*---------------------------------------------------------------------------*/
37 #include "calc.H"
39 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
41 #include "incompressible/RAS/RASModel/RASModel.H"
42 #include "incompressible/LES/LESModel/LESModel.H"
44 #include "basicPsiThermo.H"
45 #include "compressible/RAS/RASModel/RASModel.H"
46 #include "compressible/LES/LESModel/LESModel.H"
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 namespace Foam
53     void execFlowFunctionObjects(const argList& args, const Time& runTime)
54     {
55         if (args.optionFound("dict"))
56         {
57             IOdictionary dict
58             (
59                 IOobject
60                 (
61                     args.option("dict"),
62                     runTime.system(),
63                     runTime,
64                     IOobject::MUST_READ
65                 )
66             );
68             functionObjectList fol(runTime, dict);
69             fol.start();
70             fol.execute();
71         }
72         else
73         {
74             functionObjectList fol(runTime);
75             fol.start();
76             fol.execute();
77         }
78     }
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
86     Info<< "    Reading phi" << endl;
87     surfaceScalarField phi
88     (
89         IOobject
90         (
91             "phi",
92             runTime.timeName(),
93             mesh,
94             IOobject::MUST_READ
95         ),
96         mesh
97     );
99     Info<< "    Reading U" << endl;
100     volVectorField U
101     (
102         IOobject
103         (
104             "U",
105             runTime.timeName(),
106             mesh,
107             IOobject::MUST_READ
108         ),
109         mesh
110     );
112     Info<< "    Reading p" << endl;
113     volScalarField p
114     (
115         IOobject
116         (
117             "p",
118             runTime.timeName(),
119             mesh,
120             IOobject::MUST_READ
121         ),
122         mesh
123     );
125     if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
126     {
127         IOobject RASPropertiesHeader
128         (
129             "RASProperties",
130             runTime.constant(),
131             mesh,
132             IOobject::MUST_READ,
133             IOobject::NO_WRITE,
134             false
135         );
137         IOobject LESPropertiesHeader
138         (
139             "LESProperties",
140             runTime.constant(),
141             mesh,
142             IOobject::MUST_READ,
143             IOobject::NO_WRITE,
144             false
145         );
147         singlePhaseTransportModel laminarTransport(U, phi);
149         if (RASPropertiesHeader.headerOk())
150         {
151             IOdictionary RASProperties(RASPropertiesHeader);
153             autoPtr<incompressible::RASModel> RASModel
154             (
155                 incompressible::RASModel::New
156                 (
157                     U,
158                     phi,
159                     laminarTransport
160                 )
161             );
162             execFlowFunctionObjects(args, runTime);
163         }
164         else if (LESPropertiesHeader.headerOk())
165         {
166             IOdictionary LESProperties(LESPropertiesHeader);
168             autoPtr<incompressible::LESModel> sgsModel
169             (
170                 incompressible::LESModel::New(U, phi, laminarTransport)
171             );
173             execFlowFunctionObjects(args, runTime);
174         }
175         else
176         {
177             IOdictionary transportProperties
178             (
179                 IOobject
180                 (
181                     "transportProperties",
182                     runTime.constant(),
183                     mesh,
184                     IOobject::MUST_READ,
185                     IOobject::NO_WRITE
186                 )
187             );
189             dimensionedScalar nu(transportProperties.lookup("nu"));
191             execFlowFunctionObjects(args, runTime);
192         }
193     }
194     else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
195     {
196         autoPtr<basicPsiThermo> thermo(basicPsiThermo::New(mesh));
198         volScalarField rho
199         (
200             IOobject
201             (
202                 "rho",
203                 runTime.timeName(),
204                 mesh
205             ),
206             thermo->rho()
207         );
209         IOobject RASPropertiesHeader
210         (
211             "RASProperties",
212             runTime.constant(),
213             mesh,
214             IOobject::MUST_READ,
215             IOobject::NO_WRITE,
216             false
217         );
219         IOobject LESPropertiesHeader
220         (
221             "LESProperties",
222             runTime.constant(),
223             mesh,
224             IOobject::MUST_READ,
225             IOobject::NO_WRITE,
226             false
227         );
229         if (RASPropertiesHeader.headerOk())
230         {
231             IOdictionary RASProperties(RASPropertiesHeader);
233             autoPtr<compressible::RASModel> RASModel
234             (
235                 compressible::RASModel::New
236                 (
237                     rho,
238                     U,
239                     phi,
240                     thermo()
241                 )
242             );
244             execFlowFunctionObjects(args, runTime);
245         }
246         else if (LESPropertiesHeader.headerOk())
247         {
248             IOdictionary LESProperties(LESPropertiesHeader);
250             autoPtr<compressible::LESModel> sgsModel
251             (
252                 compressible::LESModel::New(rho, U, phi, thermo())
253             );
255             execFlowFunctionObjects(args, runTime);
256         }
257         else
258         {
259             IOdictionary transportProperties
260             (
261                 IOobject
262                 (
263                     "transportProperties",
264                     runTime.constant(),
265                     mesh,
266                     IOobject::MUST_READ,
267                     IOobject::NO_WRITE
268                 )
269             );
271             dimensionedScalar mu(transportProperties.lookup("mu"));
273             execFlowFunctionObjects(args, runTime);
274         }
275     }
276     else
277     {
278         FatalErrorIn(args.executable())
279             << "Incorrect dimensions of phi: " << phi.dimensions()
280             << nl << exit(FatalError);
281     }
285 // ************************************************************************* //