execFlowFunctionObjects.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  execFlowFunctionObjects
26 
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/ directory.
31 
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.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "timeSelector.H"
40 
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "pointFields.H"
45 #include "ReadFields.H"
46 #include "fvOptions.H"
47 
51 
52 using namespace Foam;
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 // Read all fields of type. Returns names of fields read. Guarantees all
57 // processors to read fields in same order.
58 template<class GeoField>
59 wordList ReadUniformFields
60 (
61  const IOobjectList& objects,
63  const bool syncPar
64 )
65 {
66  // Search list of objects for wanted type
67  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
68 
69  wordList masterNames(fieldObjects.names());
70 
71  if (syncPar && Pstream::parRun())
72  {
73  // Check that I have the same fields as the master
74  const wordList localNames(masterNames);
75  Pstream::scatter(masterNames);
76 
77  HashSet<word> localNamesSet(localNames);
78 
79  forAll(masterNames, i)
80  {
81  const word& masterFld = masterNames[i];
82 
83  HashSet<word>::iterator iter = localNamesSet.find(masterFld);
84 
85  if (iter == localNamesSet.end())
86  {
88  << "Fields not synchronised across processors." << endl
89  << "Master has fields " << masterNames
90  << " processor " << Pstream::myProcNo()
91  << " has fields " << localNames << exit(FatalError);
92  }
93  else
94  {
95  localNamesSet.erase(iter);
96  }
97  }
98 
99  forAllConstIter(HashSet<word>, localNamesSet, iter)
100  {
102  << "Fields not synchronised across processors." << endl
103  << "Master has fields " << masterNames
104  << " processor " << Pstream::myProcNo()
105  << " has fields " << localNames << exit(FatalError);
106  }
107  }
108 
109 
110  fields.setSize(masterNames.size());
111 
112  // Make sure to read in masterNames order.
113 
114  forAll(masterNames, i)
115  {
116  Info<< "Reading " << GeoField::typeName << ' ' << masterNames[i]
117  << endl;
118 
119  const IOobject& io = *fieldObjects[masterNames[i]];
120 
121  fields.set
122  (
123  i,
124  new GeoField
125  (
126  IOobject
127  (
128  io.name(),
129  io.instance(),
130  io.local(),
131  io.db(),
134  io.registerObject()
135  )
136  )
137  );
138  }
139  return masterNames;
140 }
141 
142 
143 void calc
144 (
145  const argList& args,
146  const Time& runTime,
147  const fvMesh& mesh,
148  functionObjectList& fol
149 )
150 {
151  if (args.optionFound("noRead"))
152  {
153  fol.execute(true);
154  }
155  else if (args.optionFound("noFlow"))
156  {
157  Info<< " Operating in no-flow mode; no models will be loaded."
158  << " All vol, surface and point fields will be loaded." << endl;
159 
160  // Read objects in time directory
161  IOobjectList objects(mesh, runTime.timeName());
162 
163  // Read vol fields.
164 
166  ReadFields(mesh, objects, vsFlds);
167 
169  ReadFields(mesh, objects, vvFlds);
170 
172  ReadFields(mesh, objects, vstFlds);
173 
174  PtrList<volSymmTensorField> vsymtFlds;
175  ReadFields(mesh, objects, vsymtFlds);
176 
178  ReadFields(mesh, objects, vtFlds);
179 
180  // Read vol-internal fields.
181 
183  ReadFields(mesh, objects, vsiFlds);
184 
186  ReadFields(mesh, objects, vviFlds);
187 
189  ReadFields(mesh, objects, vstiFlds);
190 
192  ReadFields(mesh, objects, vsymtiFlds);
193 
195  ReadFields(mesh, objects, vtiFlds);
196 
197  // Read surface fields.
198 
200  ReadFields(mesh, objects, ssFlds);
201 
203  ReadFields(mesh, objects, svFlds);
204 
206  ReadFields(mesh, objects, sstFlds);
207 
209  ReadFields(mesh, objects, ssymtFlds);
210 
212  ReadFields(mesh, objects, stFlds);
213 
214  // Read point fields.
215  const pointMesh& pMesh = pointMesh::New(mesh);
216 
218  ReadFields(pMesh, objects, psFlds);
219 
221  ReadFields(pMesh, objects, pvFlds);
222 
224  ReadFields(pMesh, objects, pstFlds);
225 
227  ReadFields(pMesh, objects, psymtFlds);
228 
230  ReadFields(pMesh, objects, ptFlds);
231 
232  // Read uniform dimensioned fields
233  IOobjectList constantObjects(mesh, runTime.constant());
234 
236  ReadUniformFields(constantObjects, usFlds, true);
237 
239  ReadUniformFields(constantObjects, uvFlds, true);
240 
242  ReadUniformFields(constantObjects, ustFlds, true);
243 
245  ReadUniformFields(constantObjects, usymmtFlds, true);
246 
248  ReadUniformFields(constantObjects, utFlds, true);
249 
250  fol.execute(true);
251  }
252  else
253  {
254  Info<< " Reading phi" << endl;
256  (
257  IOobject
258  (
259  "phi",
260  runTime.timeName(),
261  mesh,
263  ),
264  mesh
265  );
266 
267  Info<< " Reading U" << endl;
269  (
270  IOobject
271  (
272  "U",
273  runTime.timeName(),
274  mesh,
276  ),
277  mesh
278  );
279 
280  Info<< " Reading p" << endl;
282  (
283  IOobject
284  (
285  "p",
286  runTime.timeName(),
287  mesh,
289  ),
290  mesh
291  );
292 
293  // Note: fvOptions not directly used but constructs fvOptions so
294  // e.g. porosity modelling is effective for use in forces fo.
295  #include "createFvOptions.H"
296 
297  if (phi.dimensions() == dimVolume/dimTime)
298  {
299  IOobject turbulencePropertiesHeader
300  (
301  "turbulenceProperties",
302  runTime.constant(),
303  mesh,
306  false
307  );
308 
309  if (turbulencePropertiesHeader.headerOk())
310  {
312 
314  (
316  (
317  U,
318  phi,
320  )
321  );
322 
323  fol.execute(true);
324  }
325  else
326  {
327  IOdictionary transportProperties
328  (
329  IOobject
330  (
331  "transportProperties",
332  runTime.constant(),
333  mesh,
336  )
337  );
338 
339  fol.execute(true);
340  }
341  }
342  else if (phi.dimensions() == dimMass/dimTime)
343  {
345 
347  (
348  IOobject
349  (
350  "rho",
351  runTime.timeName(),
352  mesh
353  ),
354  thermo->rho()
355  );
356 
357  IOobject turbulencePropertiesHeader
358  (
359  "turbulenceProperties",
360  runTime.constant(),
361  mesh,
364  false
365  );
366 
367  if (turbulencePropertiesHeader.headerOk())
368  {
370  (
372  (
373  rho,
374  U,
375  phi,
376  thermo()
377  )
378  );
379 
380  fol.execute(true);
381  }
382  else
383  {
384  IOdictionary transportProperties
385  (
386  IOobject
387  (
388  "transportProperties",
389  runTime.constant(),
390  mesh,
393  )
394  );
395 
396  fol.execute(true);
397  }
398  }
399  else
400  {
402  << "Incorrect dimensions of phi: " << phi.dimensions()
403  << nl << exit(FatalError);
404  }
405  }
406 }
407 
408 
409 autoPtr<functionObjectList> readFunctionObjects
410 (
411  const argList& args,
412  const Time& runTime,
413  dictionary& folDict
414 )
415 {
417 
418  if (args.optionFound("dict"))
419  {
420  folDict = IOdictionary
421  (
422  IOobject
423  (
424  args["dict"],
425  runTime,
427  )
428  );
429  folPtr.reset(new functionObjectList(runTime, folDict));
430  }
431  else
432  {
433  folPtr.reset(new functionObjectList(runTime));
434  }
435  folPtr->start();
436 
437  return folPtr;
438 }
439 
440 
441 int main(int argc, char *argv[])
442 {
444  #include "addRegionOption.H"
446  (
447  "noFlow",
448  "suppress creating flow models"
449  );
451  (
452  "noRead",
453  "do not read any field data"
454  );
455  #include "addDictOption.H"
456 
457  #include "setRootCase.H"
458  #include "createTime.H"
460  #include "createNamedMesh.H"
461 
462  // Externally stored dictionary for functionObjectList
463  // if not constructed from runTime
464  dictionary folDict;
465 
466  // Construct functionObjectList
468  (
469  readFunctionObjects(args, runTime, folDict)
470  );
471 
472  forAll(timeDirs, timeI)
473  {
474  runTime.setTime(timeDirs[timeI], timeI);
475 
476  Info<< "Time = " << runTime.timeName() << endl;
477 
479  {
480  // Update functionObjectList if mesh changes
481  folPtr = readFunctionObjects(args, runTime, folDict);
482  }
483 
485 
486  try
487  {
488  calc(args, runTime, mesh, folPtr());
489  }
490  catch (IOerror& err)
491  {
492  Warning<< err << endl;
493  }
494 
495  Info<< endl;
496  }
497 
498  Info<< "End\n" << endl;
499 
500  return 0;
501 }
502 
503 
504 // ************************************************************************* //
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IncompressibleTurbulenceModel::New
static autoPtr< IncompressibleTurbulenceModel > New(const volVectorField &U, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
Definition: IncompressibleTurbulenceModel.C:68
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
p
p
Definition: pEqn.H:62
fvOptions.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::HashTable::iterator
An STL-conforming iterator.
Definition: HashTable.H:415
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
singlePhaseTransportModel.H
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::ReadFields
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Helper routine to read Geometric fields.
Foam::Warning
messageStream Warning
Foam::fluidThermo::New
static autoPtr< fluidThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition: fluidThermo.C:60
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::calc
void calc(const argList &args, const Time &runTime, const fvMesh &mesh)
Definition: Lambda2.C:38
turbulentTransportModel.H
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:147
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:97
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:55
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
Foam::functionObjectList
List of function objects with start(), execute() and end() functions that is called for each object.
Definition: functionObjectList.H:58
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
U
U
Definition: pEqn.H:46
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
addDictOption.H
Foam::IOobject::db
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:239
createFvOptions.H
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:498
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:60
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::IOobject::registerObject
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:303
argList.H
addRegionOption.H
Foam::polyMesh::UNCHANGED
@ UNCHANGED
Definition: polyMesh.H:90
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
createNamedMesh.H
Foam::Time::setTime
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:964
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
Foam::error::throwExceptions
void throwExceptions()
Definition: error.H:121
Foam::ThermalDiffusivity::New
static autoPtr< ThermalDiffusivity > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
Definition: ThermalDiffusivity.C:62
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
rho
rho
Definition: pEqn.H:3
uniformDimensionedFields.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
ReadFields.H
Helper routine to read fields.
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::functionObjectList::execute
virtual bool execute(const bool forceWrite=false)
Called at each ++ or += of the time-loop. forceWrite overrides.
Definition: functionObjectList.C:197
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
timeSelector.H
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::IOerror
Report an I/O error.
Definition: error.H:195
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:253
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::IOobject::local
const fileName & local() const
Definition: IOobject.H:360
pointFields.H
turbulentFluidThermoModel.H