ReadFields.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) 2011-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "ReadFields.H"
27 #include "IOobjectList.H"
28 #include "objectRegistry.H"
29 
30 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
31 
32 template<class GeoField, class Mesh>
34 (
35  const Mesh& mesh,
36  const IOobjectList& objects,
37  PtrList<GeoField>& fields,
38  const bool syncPar
39 )
40 {
41  // Search list of objects for wanted type
42  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
43 
44  wordList masterNames(fieldObjects.names());
45 
46  if (syncPar && Pstream::parRun())
47  {
48  // Check that I have the same fields as the master
49  const wordList localNames(masterNames);
50  Pstream::scatter(masterNames);
51 
52  HashSet<word> localNamesSet(localNames);
53 
54  forAll(masterNames, i)
55  {
56  const word& masterFld = masterNames[i];
57 
58  HashSet<word>::iterator iter = localNamesSet.find(masterFld);
59 
60  if (iter == localNamesSet.end())
61  {
63  << "Fields not synchronised across processors." << endl
64  << "Master has fields " << masterNames
65  << " processor " << Pstream::myProcNo()
66  << " has fields " << localNames << exit(FatalError);
67  }
68  else
69  {
70  localNamesSet.erase(iter);
71  }
72  }
73 
74  forAllConstIter(HashSet<word>, localNamesSet, iter)
75  {
77  << "Fields not synchronised across processors." << endl
78  << "Master has fields " << masterNames
79  << " processor " << Pstream::myProcNo()
80  << " has fields " << localNames << exit(FatalError);
81  }
82  }
83 
84 
85  fields.setSize(masterNames.size());
86 
87  // Make sure to read in masterNames order.
88 
89  forAll(masterNames, i)
90  {
91  Info<< "Reading " << GeoField::typeName << ' ' << masterNames[i]
92  << endl;
93 
94  const IOobject& io = *fieldObjects[masterNames[i]];
95 
96  fields.set
97  (
98  i,
99  new GeoField
100  (
101  IOobject
102  (
103  io.name(),
104  io.instance(),
105  io.local(),
106  io.db(),
107  IOobject::MUST_READ,
108  IOobject::AUTO_WRITE,
109  io.registerObject()
110  ),
111  mesh
112  )
113  );
114  }
115  return masterNames;
116 }
117 
118 
119 template<class GeoField>
120 void Foam::ReadFields
121 (
122  const word& fieldName,
123  const typename GeoField::Mesh& mesh,
124  const wordList& timeNames,
125  objectRegistry& fieldsCache
126 )
127 {
128  // Collect all times that are no longer used
129  {
130  HashSet<word> usedTimes(timeNames);
131 
132  DynamicList<word> unusedTimes(fieldsCache.size());
133 
134  forAllIter(objectRegistry, fieldsCache, timeIter)
135  {
136  const word& tm = timeIter.key();
137  if (!usedTimes.found(tm))
138  {
139  unusedTimes.append(tm);
140  }
141  }
142 
143  forAll(unusedTimes, i)
144  {
145  objectRegistry& timeCache = const_cast<objectRegistry&>
146  (
147  fieldsCache.lookupObject<objectRegistry>(unusedTimes[i])
148  );
149  fieldsCache.checkOut(timeCache);
150  }
151  }
152 
153 
154  // Load any new fields
155  forAll(timeNames, i)
156  {
157  const word& tm = timeNames[i];
158 
159  // Create if not found
160  if (!fieldsCache.found(tm))
161  {
162  // Create objectRegistry if not found
163  objectRegistry* timeCachePtr = new objectRegistry
164  (
165  IOobject
166  (
167  tm,
168  tm,
169  fieldsCache,
170  IOobject::NO_READ,
171  IOobject::NO_WRITE
172  )
173  );
174  timeCachePtr->store();
175  }
176 
177  // Obtain cache for current time
178  const objectRegistry& timeCache =
179  fieldsCache.lookupObject<objectRegistry>
180  (
181  tm
182  );
183 
184  // Store field if not found
185  if (!timeCache.found(fieldName))
186  {
187  GeoField loadedFld
188  (
189  IOobject
190  (
191  fieldName,
192  tm,
193  mesh.thisDb(),
194  IOobject::MUST_READ,
195  IOobject::NO_WRITE,
196  false
197  ),
198  mesh
199  );
200 
201  // Transfer to timeCache (new objectRegistry and store flag)
202  GeoField* fldPtr = new GeoField
203  (
204  IOobject
205  (
206  fieldName,
207  tm,
208  timeCache,
209  IOobject::NO_READ,
210  IOobject::NO_WRITE
211  ),
212  loadedFld
213  );
214  fldPtr->store();
215  }
216  }
217 }
218 
219 
220 template<class GeoField>
221 void Foam::ReadFields
222 (
223  const word& fieldName,
224  const typename GeoField::Mesh& mesh,
225  const wordList& timeNames,
226  const word& registryName
227 )
228 {
229  ReadFields<GeoField>
230  (
231  fieldName,
232  mesh,
233  timeNames,
234  const_cast<objectRegistry&>
235  (
236  mesh.thisDb().subRegistry(registryName, true)
237  )
238  );
239 }
240 
241 
242 template<class GeoFieldType>
243 void Foam::readFields
244 (
245  const typename GeoFieldType::Mesh& mesh,
246  const IOobjectList& objects,
247  const HashSet<word>& selectedFields,
248  LIFOStack<regIOobject*>& storedObjects
249 )
250 {
251  IOobjectList fields(objects.lookupClass(GeoFieldType::typeName));
252  if (!fields.size()) return;
253 
254  bool firstField = true;
255 
256  forAllConstIter(IOobjectList, fields, fieldIter)
257  {
258  const IOobject& io = *fieldIter();
259  const word& fieldName = io.name();
260 
261  if (selectedFields.found(fieldName))
262  {
263  if (firstField)
264  {
265  Info<< " " << GeoFieldType::typeName << "s:";
266  firstField = false;
267  }
268 
269  Info<< " " << fieldName;
270 
271  GeoFieldType* fieldPtr = new GeoFieldType
272  (
273  IOobject
274  (
275  fieldName,
276  io.instance(),
277  io.local(),
278  io.db(),
279  IOobject::MUST_READ,
280  IOobject::NO_WRITE
281  ),
282  mesh
283  );
284  fieldPtr->store();
285  storedObjects.push(fieldPtr);
286  }
287  }
288 
289  if (!firstField)
290  {
291  Info<< endl;
292  }
293 }
294 
295 
296 template<class UniformFieldType>
298 (
299  const IOobjectList& objects,
300  const HashSet<word>& selectedFields,
301  LIFOStack<regIOobject*>& storedObjects,
302  const bool syncPar
303 )
304 {
305  // Search list of objects for wanted type
306  IOobjectList fields(objects.lookupClass(UniformFieldType::typeName));
307  if (!fields.size()) return;
308 
309  wordList masterNames(fields.names());
310 
311  if (syncPar && Pstream::parRun())
312  {
313  // Check that I have the same fields as the master
314  const wordList localNames(masterNames);
315  Pstream::scatter(masterNames);
316 
317  HashSet<word> localNamesSet(localNames);
318 
319  forAll(masterNames, i)
320  {
321  const word& masterFld = masterNames[i];
322 
323  HashSet<word>::iterator iter = localNamesSet.find(masterFld);
324 
325  if (iter == localNamesSet.end())
326  {
328  << "Fields not synchronised across processors." << endl
329  << "Master has fields " << masterNames
330  << " processor " << Pstream::myProcNo()
331  << " has fields " << localNames << exit(FatalError);
332  }
333  else
334  {
335  localNamesSet.erase(iter);
336  }
337  }
338 
339  forAllConstIter(HashSet<word>, localNamesSet, iter)
340  {
342  << "Fields not synchronised across processors." << endl
343  << "Master has fields " << masterNames
344  << " processor " << Pstream::myProcNo()
345  << " has fields " << localNames << exit(FatalError);
346  }
347  }
348 
349  bool firstField = true;
350 
351  forAll(masterNames, i)
352  {
353  const IOobject& io = *fields[masterNames[i]];
354  const word& fieldName = io.name();
355 
356  if (selectedFields.found(fieldName))
357  {
358  if (firstField)
359  {
360  Info<< " " << UniformFieldType::typeName << "s:";
361  firstField = false;
362  }
363 
364  Info<< " " << fieldName;
365 
366  UniformFieldType* fieldPtr = new UniformFieldType
367  (
368  IOobject
369  (
370  fieldName,
371  io.instance(),
372  io.local(),
373  io.db(),
374  IOobject::MUST_READ,
375  IOobject::NO_WRITE
376  )
377  );
378  fieldPtr->store();
379  storedObjects.push(fieldPtr);
380  }
381  }
382 
383  Info<< endl;
384 }
385 
386 
387 // ************************************************************************* //
A LIFO stack based on a singly-linked list.
Definition: LIFOStack.H:51
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
const word & name() const
Return name.
Definition: IOobject.H:295
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Field reading functions for post-processing utilities.
void push(const T &a)
Push an element onto the stack.
Definition: LIFOStack.H:84
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
An STL-conforming iterator.
Definition: HashTable.H:426
dynamicFvMesh & mesh
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
const fileName & local() const
Definition: IOobject.H:388
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
List< word > wordList
A List of words.
Definition: fileName.H:54
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:192
const fileName & instance() const
Definition: IOobject.H:378
messageStream Info
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:354
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void readUniformFields(const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects, const bool syncPar=true)
Read the selected UniformDimensionedFields of the specified type.
Definition: ReadFields.C:298