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-------------------------------------------------------------------------------
11License
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
35template<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
92template<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
145template<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
196template<class GeoField>
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 )
241 );
242 timeCachePtr->store();
243 }
244
245 // Obtain cache for current time
246 const objectRegistry& timeCache =
247 fieldsCache.lookupObject<objectRegistry>(timeName);
248
249 // Store field if not found
250 if (!timeCache.found(fieldName))
251 {
252 //Info<< "Loading field " << fieldName
253 // << " for time " << timeName << endl;
254
255 GeoField loadedFld
256 (
257 IOobject
258 (
259 fieldName,
260 timeName,
261 mesh.thisDb(),
262 IOobject::MUST_READ,
263 IOobject::NO_WRITE,
264 false // do not register
265 ),
266 mesh
267 );
268
269 // Transfer to timeCache (new objectRegistry and store flag)
270 GeoField* fldPtr = new GeoField
271 (
272 IOobject
273 (
274 fieldName,
275 timeName,
276 timeCache,
277 IOobject::NO_READ,
278 IOobject::NO_WRITE
279 ),
280 loadedFld
281 );
282 fldPtr->store();
283 }
284 }
285}
286
287
288template<class GeoField>
290(
291 const word& fieldName,
292 const typename GeoField::Mesh& mesh,
293 const wordList& timeNames,
294 const word& registryName
295)
296{
297 ReadFields<GeoField>
298 (
299 fieldName,
300 mesh,
301 timeNames,
302 const_cast<objectRegistry&>
303 (
304 mesh.thisDb().subRegistry(registryName, true)
305 )
306 );
307}
308
309
310template<class GeoFieldType>
312(
313 const typename GeoFieldType::Mesh& mesh,
314 const IOobjectList& objects,
315 const wordHashSet& selectedFields,
316 LIFOStack<regIOobject*>& storedObjects
317)
318{
319 // Names of GeoField objects, sorted order. Not synchronised.
320 const wordList fieldNames
321 (
322 objects.sortedNames
323 (
324 GeoFieldType::typeName,
325 selectedFields // Only permit these
326 )
327 );
328
329 label nFields = 0;
330
331 for (const word& fieldName : fieldNames)
332 {
333 const IOobject& io = *objects[fieldName];
334
335 if (!nFields)
336 {
337 Info<< " " << GeoFieldType::typeName << ':';
338 }
339 Info<< ' ' << fieldName;
340
341 GeoFieldType* fieldPtr = new GeoFieldType
342 (
344 (
345 fieldName,
346 io.instance(),
347 io.local(),
348 io.db(),
349 IOobject::MUST_READ,
350 IOobject::NO_WRITE
351 ),
352 mesh
353 );
354 fieldPtr->store();
355 storedObjects.push(fieldPtr);
356
357 ++nFields;
358 }
359
360 if (nFields) Info<< endl;
361}
362
363
364template<class UniformFieldType>
366(
367 const IOobjectList& objects,
368 const wordHashSet& selectedFields,
369 LIFOStack<regIOobject*>& storedObjects,
370 const bool syncPar
371)
372{
373 // Names of UniformField objects, sorted order.
374 const wordList fieldNames
375 (
376 objects.names
377 (
378 UniformFieldType::typeName,
379 selectedFields, // Only permit these
380 syncPar
381 )
382 );
383
384 label nFields = 0;
385
386 for (const word& fieldName : fieldNames)
387 {
388 const IOobject& io = *objects[fieldName];
389
390 if (!nFields)
391 {
392 Info<< " " << UniformFieldType::typeName << ':';
393 }
394 Info<< ' ' << fieldName;
395
396 UniformFieldType* fieldPtr = new UniformFieldType
397 (
399 (
400 fieldName,
401 io.instance(),
402 io.local(),
403 io.db(),
404 IOobject::MUST_READ,
405 IOobject::NO_WRITE
406 )
407 );
408 fieldPtr->store();
409 storedObjects.push(fieldPtr);
410
411 ++nFields;
412 }
413
414 if (nFields) Info<< endl;
415}
416
417
418// ************************************************************************* //
Field reading functions for post-processing utilities.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:383
wordList names() const
The unsorted names of the IOobjects.
Definition: IOobjectList.C:351
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A LIFO stack based on a singly-linked list.
Definition: LIFOStack.H:54
void push(const T &elem)
Push an element onto the front of the stack.
Definition: LIFOStack.H:78
faMesh Mesh
The mesh type.
Definition: faMesh.H:491
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
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.
List< word > wordList
A List of words.
Definition: fileName.H:63
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
void readUniformFields(const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects, const bool syncPar=true)
Read the selected UniformDimensionedFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:82
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
static const char *const typeName
The type name used in ensight case files.