ReadFieldsTemplates.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2014 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "ReadFields.H"
30 #include "HashSet.H"
31 #include "IOobjectList.H"
32 
33 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
34 
35 template<class Type, template<class> class PatchField, class GeoMesh>
37 (
38  const typename GeoMesh::Mesh& mesh,
39  const IOobjectList& objects,
40  PtrList<GeometricField<Type, PatchField, GeoMesh>>& fields,
41  const bool syncPar,
42  const bool readOldTime
43 )
44 {
45  typedef GeometricField<Type, PatchField, GeoMesh> GeoField;
46 
47  // Names of GeoField objects, sorted order.
48  const wordList fieldNames(objects.names(GeoField::typeName, syncPar));
49 
50  // Construct the fields - reading in consistent (master) order.
51  fields.resize(fieldNames.size());
52 
53  label nFields = 0;
54 
55  for (const word& fieldName : fieldNames)
56  {
57  if (!nFields)
58  {
59  Info<< "Reading " << GeoField::typeName << ':';
60  }
61  Info<< ' ' << fieldName;
62 
63  const IOobject& io = *objects[fieldName];
64 
65  fields.set
66  (
67  nFields++,
68  new GeoField
69  (
70  IOobject
71  (
72  io.name(),
73  io.instance(),
74  io.local(),
75  io.db(),
76  IOobject::MUST_READ,
77  IOobject::AUTO_WRITE,
78  io.registerObject()
79  ),
80  mesh,
81  readOldTime
82  )
83  );
84  }
85 
86  if (nFields) Info<< endl;
87 
88  return fieldNames;
89 }
90 
91 
92 template<class GeoField, class Mesh>
94 (
95  const Mesh& mesh,
96  const IOobjectList& objects,
97  PtrList<GeoField>& fields,
98  const bool syncPar
99 )
100 {
101  // Names of GeoField objects, sorted order.
102  const wordList fieldNames(objects.names(GeoField::typeName, syncPar));
103 
104  // Construct the fields - reading in consistent (master) order.
105  fields.resize(fieldNames.size());
106 
107  label nFields = 0;
108 
109  for (const word& fieldName : fieldNames)
110  {
111  if (!nFields)
112  {
113  Info<< "Reading " << GeoField::typeName << ':';
114  }
115  Info<< ' ' << fieldName;
116 
117  const IOobject& io = *objects[fieldName];
118 
119  fields.set
120  (
121  nFields++,
122  new GeoField
123  (
124  IOobject
125  (
126  io.name(),
127  io.instance(),
128  io.local(),
129  io.db(),
130  IOobject::MUST_READ,
131  IOobject::AUTO_WRITE,
132  io.registerObject()
133  ),
134  mesh
135  )
136  );
137  }
138 
139  if (nFields) Info<< endl;
140 
141  return fieldNames;
142 }
143 
144 
145 template<class GeoField>
147 (
148  const IOobjectList& objects,
149  PtrList<GeoField>& fields,
150  const bool syncPar
151 )
152 {
153  // Names of GeoField objects, sorted order.
154  const wordList fieldNames(objects.names(GeoField::typeName, syncPar));
155 
156  // Construct the fields - reading in consistent (master) order.
157  fields.resize(fieldNames.size());
158 
159  label nFields = 0;
160 
161  for (const word& fieldName : fieldNames)
162  {
163  if (!nFields)
164  {
165  Info<< "Reading " << GeoField::typeName << ':';
166  }
167  Info<< ' ' << fieldName;
168 
169  const IOobject& io = *objects[fieldName];
170 
171  fields.set
172  (
173  nFields++,
174  new GeoField
175  (
176  IOobject
177  (
178  io.name(),
179  io.instance(),
180  io.local(),
181  io.db(),
182  IOobject::MUST_READ,
183  IOobject::AUTO_WRITE,
184  io.registerObject()
185  )
186  )
187  );
188  }
189 
190  if (nFields) Info<< endl;
191 
192  return fieldNames;
193 }
194 
195 
196 template<class GeoField>
197 void Foam::ReadFields
198 (
199  const word& fieldName,
200  const typename GeoField::Mesh& mesh,
201  const wordList& timeNames,
202  objectRegistry& fieldsCache
203 )
204 {
205  // Unload times that are no longer used
206  {
207  wordHashSet unusedTimes(fieldsCache.toc());
208  unusedTimes.erase(timeNames);
209 
210  //Info<< "Unloading times " << unusedTimes << endl;
211 
212  for (const word& timeName : unusedTimes)
213  {
214  objectRegistry& timeCache =
215  fieldsCache.lookupObjectRef<objectRegistry>(timeName);
216 
217  fieldsCache.checkOut(timeCache);
218  }
219  }
220 
221 
222  // Load any new fields
223  for (const word& timeName : timeNames)
224  {
225  // Create if not found
226  if (!fieldsCache.found(timeName))
227  {
228  //Info<< "Creating registry for time " << timeName << endl;
229 
230  // Create objectRegistry if not found
231  objectRegistry* timeCachePtr = new objectRegistry
232  (
233  IOobject
234  (
235  timeName,
236  timeName,
237  fieldsCache,
238  IOobject::NO_READ,
239  IOobject::NO_WRITE,
240  IOobject::REGISTER
241  )
242  );
243  timeCachePtr->store();
244  }
245 
246  // Obtain cache for current time
247  const objectRegistry& timeCache =
248  fieldsCache.lookupObject<objectRegistry>(timeName);
249 
250  // Store field if not found
251  if (!timeCache.found(fieldName))
252  {
253  //Info<< "Loading field " << fieldName
254  // << " for time " << timeName << endl;
255 
256  GeoField loadedFld
257  (
258  IOobject
259  (
260  fieldName,
261  timeName,
262  mesh.thisDb(),
263  IOobject::MUST_READ,
264  IOobject::NO_WRITE,
265  IOobject::NO_REGISTER
266  ),
267  mesh
268  );
269 
270  // Transfer to timeCache (new objectRegistry and store flag)
271  GeoField* fldPtr = new GeoField
272  (
273  IOobject
274  (
275  fieldName,
276  timeName,
277  timeCache,
278  IOobject::NO_READ,
279  IOobject::NO_WRITE,
280  IOobject::REGISTER
281  ),
282  loadedFld
283  );
284  fldPtr->store();
285  }
286  }
287 }
288 
289 
290 template<class GeoField>
291 void Foam::ReadFields
292 (
293  const word& fieldName,
294  const typename GeoField::Mesh& mesh,
295  const wordList& timeNames,
296  const word& registryName
297 )
298 {
299  ReadFields<GeoField>
300  (
301  fieldName,
302  mesh,
303  timeNames,
304  const_cast<objectRegistry&>
305  (
306  mesh.thisDb().subRegistry(registryName, true)
307  )
308  );
309 }
310 
311 
312 template<class GeoFieldType>
313 void Foam::readFields
314 (
315  const typename GeoFieldType::Mesh& mesh,
316  const IOobjectList& objects,
317  const wordHashSet& selectedFields,
318  LIFOStack<regIOobject*>& storedObjects
319 )
320 {
321  // Names of GeoField objects, sorted order. Not synchronised.
322  const wordList fieldNames
323  (
324  objects.sortedNames
325  (
326  GeoFieldType::typeName,
327  selectedFields // Only permit these
328  )
329  );
330 
331  label nFields = 0;
332 
333  for (const word& fieldName : fieldNames)
334  {
335  const IOobject& io = *objects[fieldName];
336 
337  if (!nFields)
338  {
339  Info<< " " << GeoFieldType::typeName << ':';
340  }
341  Info<< ' ' << fieldName;
342 
343  GeoFieldType* fieldPtr = new GeoFieldType
344  (
345  IOobject
346  (
347  fieldName,
348  io.instance(),
349  io.local(),
350  io.db(),
351  IOobject::MUST_READ,
352  IOobject::NO_WRITE,
353  IOobject::REGISTER
354  ),
355  mesh
356  );
357  fieldPtr->store();
358  storedObjects.push(fieldPtr);
359 
360  ++nFields;
361  }
363  if (nFields) Info<< endl;
364 }
365 
366 
367 template<class UniformFieldType>
369 (
370  const IOobjectList& objects,
371  const wordHashSet& selectedFields,
372  LIFOStack<regIOobject*>& storedObjects,
373  const bool syncPar
374 )
375 {
376  // Names of UniformField objects, sorted order.
377  const wordList fieldNames
378  (
379  objects.names
380  (
381  UniformFieldType::typeName,
382  selectedFields, // Only permit these
383  syncPar
384  )
385  );
386 
387  label nFields = 0;
388 
389  for (const word& fieldName : fieldNames)
390  {
391  const IOobject& io = *objects[fieldName];
392 
393  if (!nFields)
394  {
395  Info<< " " << UniformFieldType::typeName << ':';
396  }
397  Info<< ' ' << fieldName;
398 
399  UniformFieldType* fieldPtr = new UniformFieldType
400  (
401  IOobject
402  (
403  fieldName,
404  io.instance(),
405  io.local(),
406  io.db(),
407  IOobject::MUST_READ,
408  IOobject::NO_WRITE,
409  IOobject::REGISTER
410  )
411  );
412  fieldPtr->store();
413  storedObjects.push(fieldPtr);
414 
415  ++nFields;
416  }
417 
418  if (nFields) Info<< endl;
419 }
420 
421 
422 // ************************************************************************* //
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Field reading functions for post-processing utilities.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:473
List< word > wordList
List of word.
Definition: fileName.H:58
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readUniformFields(const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects, const bool syncPar=true)
Read the selected UniformDimensionedFields of the templated type.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)