1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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.
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 \*---------------------------------------------------------------------------*/
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 void execFlowFunctionObjects(const argList& args, const Time& runTime)
55 if (args.optionFound("dict"))
68 functionObjectList fol(runTime, dict);
74 functionObjectList fol(runTime);
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
86 Info<< " Reading phi" << endl;
87 surfaceScalarField phi
99 Info<< " Reading U" << endl;
112 Info<< " Reading p" << endl;
125 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
127 IOobject RASPropertiesHeader
137 IOobject LESPropertiesHeader
147 singlePhaseTransportModel laminarTransport(U, phi);
149 if (RASPropertiesHeader.headerOk())
151 IOdictionary RASProperties(RASPropertiesHeader);
153 autoPtr<incompressible::RASModel> RASModel
155 incompressible::RASModel::New
162 execFlowFunctionObjects(args, runTime);
164 else if (LESPropertiesHeader.headerOk())
166 IOdictionary LESProperties(LESPropertiesHeader);
168 autoPtr<incompressible::LESModel> sgsModel
170 incompressible::LESModel::New(U, phi, laminarTransport)
173 execFlowFunctionObjects(args, runTime);
177 IOdictionary transportProperties
181 "transportProperties",
189 dimensionedScalar nu(transportProperties.lookup("nu"));
191 execFlowFunctionObjects(args, runTime);
194 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
196 autoPtr<basicPsiThermo> thermo(basicPsiThermo::New(mesh));
209 IOobject RASPropertiesHeader
219 IOobject LESPropertiesHeader
229 if (RASPropertiesHeader.headerOk())
231 IOdictionary RASProperties(RASPropertiesHeader);
233 autoPtr<compressible::RASModel> RASModel
235 compressible::RASModel::New
244 execFlowFunctionObjects(args, runTime);
246 else if (LESPropertiesHeader.headerOk())
248 IOdictionary LESProperties(LESPropertiesHeader);
250 autoPtr<compressible::LESModel> sgsModel
252 compressible::LESModel::New(rho, U, phi, thermo())
255 execFlowFunctionObjects(args, runTime);
259 IOdictionary transportProperties
263 "transportProperties",
271 dimensionedScalar mu(transportProperties.lookup("mu"));
273 execFlowFunctionObjects(args, runTime);
278 FatalErrorIn(args.executable())
279 << "Incorrect dimensions of phi: " << phi.dimensions()
280 << nl << exit(FatalError);
285 // ************************************************************************* //