ENH: execFlowFunctionObjects: added -region option. Added -noFlow option.
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / miscellaneous / execFlowFunctionObjects / execFlowFunctionObjects.C
blobf0487ac2a2a5c5c9b62ff39ba8168c93a1615099
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 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.
30     Alternative dictionaries should be placed in the system/ folder.
32     The flow (p-U) and optionally turbulence fields are available for the
33     function objects to operate on allowing forces and other related properties
34     to be calculated in addition to cutting planes etc.
36 \*---------------------------------------------------------------------------*/
38 #include "calc.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "pointFields.H"
43 #include "ReadFields.H"
45 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
47 #include "incompressible/RAS/RASModel/RASModel.H"
48 #include "incompressible/LES/LESModel/LESModel.H"
50 #include "basicPsiThermo.H"
51 #include "compressible/RAS/RASModel/RASModel.H"
52 #include "compressible/LES/LESModel/LESModel.H"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 namespace Foam
59     void execFlowFunctionObjects(const argList& args, const Time& runTime)
60     {
61         if (args.optionFound("dict"))
62         {
63             IOdictionary dict
64             (
65                 IOobject
66                 (
67                     args["dict"],
68                     runTime.system(),
69                     runTime,
70                     IOobject::MUST_READ_IF_MODIFIED
71                 )
72             );
74             functionObjectList fol(runTime, dict);
75             fol.start();
76             fol.execute(true);  // override outputControl - force writing
77         }
78         else
79         {
80             functionObjectList fol(runTime);
81             fol.start();
82             fol.execute(true);  // override outputControl - force writing
83         }
84     }
88 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
92     if (args.optionFound("noFlow"))
93     {
94         Info<< "    Operating in no-flow mode; no models will be loaded."
95             << " All vol, surface and point fields will be loaded." << endl;
97         // Read objects in time directory
98         IOobjectList objects(mesh, runTime.timeName());
100         // Read vol fields.
102         PtrList<volScalarField> vsFlds;
103         ReadFields(mesh, objects, vsFlds);
105         PtrList<volVectorField> vvFlds;
106         ReadFields(mesh, objects, vvFlds);
108         PtrList<volSphericalTensorField> vstFlds;
109         ReadFields(mesh, objects, vstFlds);
111         PtrList<volSymmTensorField> vsymtFlds;
112         ReadFields(mesh, objects, vsymtFlds);
114         PtrList<volTensorField> vtFlds;
115         ReadFields(mesh, objects, vtFlds);
117         // Read surface fields.
119         PtrList<surfaceScalarField> ssFlds;
120         ReadFields(mesh, objects, ssFlds);
122         PtrList<surfaceVectorField> svFlds;
123         ReadFields(mesh, objects, svFlds);
125         PtrList<surfaceSphericalTensorField> sstFlds;
126         ReadFields(mesh, objects, sstFlds);
128         PtrList<surfaceSymmTensorField> ssymtFlds;
129         ReadFields(mesh, objects, ssymtFlds);
131         PtrList<surfaceTensorField> stFlds;
132         ReadFields(mesh, objects, stFlds);
134         // Read point fields.
135         const pointMesh& pMesh = pointMesh::New(mesh);
137         PtrList<pointScalarField> psFlds;
138         ReadFields(pMesh, objects, psFlds);
140         PtrList<pointVectorField> pvFlds;
141         ReadFields(pMesh, objects, pvFlds);
143         PtrList<pointSphericalTensorField> pstFlds;
144         ReadFields(pMesh, objects, pstFlds);
146         PtrList<pointSymmTensorField> psymtFlds;
147         ReadFields(pMesh, objects, psymtFlds);
149         PtrList<pointTensorField> ptFlds;
150         ReadFields(pMesh, objects, ptFlds);
152         execFlowFunctionObjects(args, runTime);
153     }
154     else
155     {
156         Info<< "    Reading phi" << endl;
157         surfaceScalarField phi
158         (
159             IOobject
160             (
161                 "phi",
162                 runTime.timeName(),
163                 mesh,
164                 IOobject::MUST_READ
165             ),
166             mesh
167         );
169         Info<< "    Reading U" << endl;
170         volVectorField U
171         (
172             IOobject
173             (
174                 "U",
175                 runTime.timeName(),
176                 mesh,
177                 IOobject::MUST_READ
178             ),
179             mesh
180         );
182         Info<< "    Reading p" << endl;
183         volScalarField p
184         (
185             IOobject
186             (
187                 "p",
188                 runTime.timeName(),
189                 mesh,
190                 IOobject::MUST_READ
191             ),
192             mesh
193         );
195         if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
196         {
197             IOobject RASPropertiesHeader
198             (
199                 "RASProperties",
200                 runTime.constant(),
201                 mesh,
202                 IOobject::MUST_READ_IF_MODIFIED,
203                 IOobject::NO_WRITE,
204                 false
205             );
207             IOobject LESPropertiesHeader
208             (
209                 "LESProperties",
210                 runTime.constant(),
211                 mesh,
212                 IOobject::MUST_READ_IF_MODIFIED,
213                 IOobject::NO_WRITE,
214                 false
215             );
217             if (RASPropertiesHeader.headerOk())
218             {
219                 IOdictionary RASProperties(RASPropertiesHeader);
221                 singlePhaseTransportModel laminarTransport(U, phi);
223                 autoPtr<incompressible::RASModel> RASModel
224                 (
225                     incompressible::RASModel::New
226                     (
227                         U,
228                         phi,
229                         laminarTransport
230                     )
231                 );
232                 execFlowFunctionObjects(args, runTime);
233             }
234             else if (LESPropertiesHeader.headerOk())
235             {
236                 IOdictionary LESProperties(LESPropertiesHeader);
238                 singlePhaseTransportModel laminarTransport(U, phi);
240                 autoPtr<incompressible::LESModel> sgsModel
241                 (
242                     incompressible::LESModel::New(U, phi, laminarTransport)
243                 );
245                 execFlowFunctionObjects(args, runTime);
246             }
247             else
248             {
249                 IOdictionary transportProperties
250                 (
251                     IOobject
252                     (
253                         "transportProperties",
254                         runTime.constant(),
255                         mesh,
256                         IOobject::MUST_READ_IF_MODIFIED,
257                         IOobject::NO_WRITE
258                     )
259                 );
261                 dimensionedScalar nu(transportProperties.lookup("nu"));
263                 execFlowFunctionObjects(args, runTime);
264             }
265         }
266         else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
267         {
268             autoPtr<basicPsiThermo> thermo(basicPsiThermo::New(mesh));
270             volScalarField rho
271             (
272                 IOobject
273                 (
274                     "rho",
275                     runTime.timeName(),
276                     mesh
277                 ),
278                 thermo->rho()
279             );
281             IOobject RASPropertiesHeader
282             (
283                 "RASProperties",
284                 runTime.constant(),
285                 mesh,
286                 IOobject::MUST_READ_IF_MODIFIED,
287                 IOobject::NO_WRITE,
288                 false
289             );
291             IOobject LESPropertiesHeader
292             (
293                 "LESProperties",
294                 runTime.constant(),
295                 mesh,
296                 IOobject::MUST_READ_IF_MODIFIED,
297                 IOobject::NO_WRITE,
298                 false
299             );
301             if (RASPropertiesHeader.headerOk())
302             {
303                 IOdictionary RASProperties(RASPropertiesHeader);
305                 autoPtr<compressible::RASModel> RASModel
306                 (
307                     compressible::RASModel::New
308                     (
309                         rho,
310                         U,
311                         phi,
312                         thermo()
313                     )
314                 );
316                 execFlowFunctionObjects(args, runTime);
317             }
318             else if (LESPropertiesHeader.headerOk())
319             {
320                 IOdictionary LESProperties(LESPropertiesHeader);
322                 autoPtr<compressible::LESModel> sgsModel
323                 (
324                     compressible::LESModel::New(rho, U, phi, thermo())
325                 );
327                 execFlowFunctionObjects(args, runTime);
328             }
329             else
330             {
331                 IOdictionary transportProperties
332                 (
333                     IOobject
334                     (
335                         "transportProperties",
336                         runTime.constant(),
337                         mesh,
338                         IOobject::MUST_READ_IF_MODIFIED,
339                         IOobject::NO_WRITE
340                     )
341                 );
343                 dimensionedScalar mu(transportProperties.lookup("mu"));
345                 execFlowFunctionObjects(args, runTime);
346             }
347         }
348         else
349         {
350             FatalErrorIn(args.executable())
351                 << "Incorrect dimensions of phi: " << phi.dimensions()
352                 << nl << exit(FatalError);
353         }
354     }
358 // ************************************************************************* //