postProcess.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2022 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  postProcess
26 
27 Description
28  Execute the set of functionObjects specified in the selected dictionary
29  (which defaults to system/controlDict) or on the command-line for the
30  selected set of times on the selected set of fields.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "argList.H"
35 #include "timeSelector.H"
36 #include "ReadFields.H"
37 #include "volFields.H"
38 #include "surfaceFields.H"
39 #include "pointFields.H"
41 
42 using namespace Foam;
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 #define ReadFields(GeoFieldType) \
47  readFields<GeoFieldType>(mesh, objects, requiredFields, storedObjects);
48 
49 #define ReadPointFields(GeoFieldType) \
50  readFields<GeoFieldType>(pMesh, objects, requiredFields, storedObjects);
51 
52 #define ReadUniformFields(FieldType) \
53  readUniformFields<FieldType> \
54  (constantObjects, requiredFields, storedObjects);
55 
56 void executeFunctionObjects
57 (
58  const argList& args,
59  const Time& runTime,
60  fvMesh& mesh,
61  const HashSet<word>& requiredFields0,
62  functionObjectList& functions,
63  bool lastTime
64 )
65 {
66  Info<< nl << "Reading fields:" << endl;
67 
68  // Maintain a stack of the stored objects to clear after executing
69  // the functionObjects
70  LIFOStack<regIOobject*> storedObjects;
71 
72  // Read objects in time directory
73  IOobjectList objects(mesh, runTime.timeName());
74 
75  HashSet<word> requiredFields(requiredFields0);
76  forAll(functions, i)
77  {
78  requiredFields.insert(functions[i].fields());
79  }
80 
81  // Read volFields
87 
88  // Read internal fields
94 
95  // Read surface fields
101 
102  // Read point fields.
103  const pointMesh& pMesh = pointMesh::New(mesh);
104 
105  ReadPointFields(pointScalarField)
106  ReadPointFields(pointVectorField);
107  ReadPointFields(pointSphericalTensorField);
108  ReadPointFields(pointSymmTensorField);
109  ReadPointFields(pointTensorField);
110 
111  // Read uniform dimensioned fields
112  IOobjectList constantObjects(mesh, runTime.constant());
113 
114  ReadUniformFields(uniformDimensionedScalarField);
115  ReadUniformFields(uniformDimensionedVectorField);
116  ReadUniformFields(uniformDimensionedSphericalTensorField);
117  ReadUniformFields(uniformDimensionedSymmTensorField);
118  ReadUniformFields(uniformDimensionedTensorField);
119 
120  Info<< nl << "Executing functionObjects" << endl;
121 
122  // Execute the functionObjects in post-processing mode
123  functions.execute();
124 
125  // Execute the functionObject 'end()' function for the last time
126  if (lastTime)
127  {
128  functions.end();
129  }
130 
131  while (!storedObjects.empty())
132  {
133  storedObjects.pop()->checkOut();
134  }
135 }
136 
137 
138 int main(int argc, char *argv[])
139 {
141  #include "addRegionOption.H"
142  #include "addFunctionObjectOptions.H"
143 
144  // Set functionObject post-processing mode
146 
147  #include "setRootCase.H"
148 
149  if (args.optionFound("list"))
150  {
151  Info<< nl
152  << "Available configured functionObjects:"
154  << endl;
155  return 0;
156  }
157 
158  #include "createTime.H"
159  Foam::instantList timeDirs = Foam::timeSelector::select0(runTime, args);
160  #include "createNamedMesh.H"
161 
162  // Initialise the set of selected fields from the command-line options
163  HashSet<word> requiredFields;
164  if (args.optionFound("fields"))
165  {
166  args.optionLookup("fields")() >> requiredFields;
167  }
168  if (args.optionFound("field"))
169  {
170  requiredFields.insert(args.optionLookup("field")());
171  }
172 
173  // Externally stored dictionary for functionObjectList
174  // if not constructed from runTime
175  dictionary functionsControlDict("controlDict");
176 
177  // Construct functionObjectList
178  autoPtr<functionObjectList> functionsPtr
179  (
181  (
182  args,
183  runTime,
184  functionsControlDict
185  )
186  );
187 
188  forAll(timeDirs, timei)
189  {
190  runTime.setTime(timeDirs[timei], timei);
191 
192  Info<< "Time = " << runTime.userTimeName() << endl;
193 
194  if (mesh.readUpdate() != polyMesh::UNCHANGED)
195  {
196  // Update functionObjectList if mesh changes
197  functionsPtr = functionObjectList::New
198  (
199  args,
200  runTime,
201  functionsControlDict
202  );
203  }
204 
206 
207  try
208  {
209  executeFunctionObjects
210  (
211  args,
212  runTime,
213  mesh,
214  requiredFields,
215  functionsPtr(),
216  timei == timeDirs.size()-1
217  );
218  }
219  catch (IOerror& err)
220  {
221  Warning<< err << endl;
222  }
223 
224  Info<< endl;
225  }
226 
227  Info<< "End\n" << endl;
228 
229  return 0;
230 }
231 
232 
233 // ************************************************************************* //
A LIFO stack based on a singly-linked list.
Definition: LIFOStack.H:51
Foam::surfaceFields.
A HashTable with keys but without contents.
Definition: HashSet.H:59
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
static wordList list()
Return the list of functionObject configuration files in.
static bool postProcess
Global post-processing mode switch.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
Field reading functions for post-processing utilities.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
bool execute()
Called at each ++ or += of the time-loop.
T pop()
Pop the top element off the stack.
Definition: LIFOStack.H:90
void throwExceptions()
Definition: error.H:115
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:959
List of function objects with start(), execute() and end() functions that is called for each object...
Report an I/O error.
Definition: error.H:187
bool checkOut()
Remove object from registry.
Definition: regIOobject.C:241
static const char nl
Definition: Ostream.H:260
objects
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
readUpdateState readUpdate(const stitchType stitch=stitchType::geometric)
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:713
messageStream Warning
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
static autoPtr< functionObjectList > New(const argList &args, const Time &runTime, dictionary &controlDict)
Construct and return a functionObjectList for an application.
word userTimeName() const
Return current user time name with units.
Definition: Time.C:863
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
bool end()
Called when Time::run() determines that the time-loop exits.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Namespace for OpenFOAM.
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
IOerror FatalIOError