1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
25 execFlowFunctionObjects
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 \*---------------------------------------------------------------------------*/
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 void execFlowFunctionObjects(const argList& args, const Time& runTime)
61 if (args.optionFound("dict"))
70 IOobject::MUST_READ_IF_MODIFIED
74 functionObjectList fol(runTime, dict);
76 fol.execute(true); // override outputControl - force writing
80 functionObjectList fol(runTime);
82 fol.execute(true); // override outputControl - force writing
88 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
92 if (args.optionFound("noFlow"))
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());
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);
156 Info<< " Reading phi" << endl;
157 surfaceScalarField phi
169 Info<< " Reading U" << endl;
182 Info<< " Reading p" << endl;
195 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
197 IOobject RASPropertiesHeader
202 IOobject::MUST_READ_IF_MODIFIED,
207 IOobject LESPropertiesHeader
212 IOobject::MUST_READ_IF_MODIFIED,
217 if (RASPropertiesHeader.headerOk())
219 IOdictionary RASProperties(RASPropertiesHeader);
221 singlePhaseTransportModel laminarTransport(U, phi);
223 autoPtr<incompressible::RASModel> RASModel
225 incompressible::RASModel::New
232 execFlowFunctionObjects(args, runTime);
234 else if (LESPropertiesHeader.headerOk())
236 IOdictionary LESProperties(LESPropertiesHeader);
238 singlePhaseTransportModel laminarTransport(U, phi);
240 autoPtr<incompressible::LESModel> sgsModel
242 incompressible::LESModel::New(U, phi, laminarTransport)
245 execFlowFunctionObjects(args, runTime);
249 IOdictionary transportProperties
253 "transportProperties",
256 IOobject::MUST_READ_IF_MODIFIED,
261 dimensionedScalar nu(transportProperties.lookup("nu"));
263 execFlowFunctionObjects(args, runTime);
266 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
268 autoPtr<basicPsiThermo> thermo(basicPsiThermo::New(mesh));
281 IOobject RASPropertiesHeader
286 IOobject::MUST_READ_IF_MODIFIED,
291 IOobject LESPropertiesHeader
296 IOobject::MUST_READ_IF_MODIFIED,
301 if (RASPropertiesHeader.headerOk())
303 IOdictionary RASProperties(RASPropertiesHeader);
305 autoPtr<compressible::RASModel> RASModel
307 compressible::RASModel::New
316 execFlowFunctionObjects(args, runTime);
318 else if (LESPropertiesHeader.headerOk())
320 IOdictionary LESProperties(LESPropertiesHeader);
322 autoPtr<compressible::LESModel> sgsModel
324 compressible::LESModel::New(rho, U, phi, thermo())
327 execFlowFunctionObjects(args, runTime);
331 IOdictionary transportProperties
335 "transportProperties",
338 IOobject::MUST_READ_IF_MODIFIED,
343 dimensionedScalar mu(transportProperties.lookup("mu"));
345 execFlowFunctionObjects(args, runTime);
350 FatalErrorIn(args.executable())
351 << "Incorrect dimensions of phi: " << phi.dimensions()
352 << nl << exit(FatalError);
358 // ************************************************************************* //