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 | Copyright (C) 2011-2014 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 "HashSet.H"
28 #include "IOobjectList.H"
29 
30 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
31 
32 // Read all GeometricFields of type. Returns names of fields read. Guarantees
33 // all processors to read fields in same order.
34 template<class Type, template<class> class PatchField, class GeoMesh>
36 (
37  const typename GeoMesh::Mesh& mesh,
38  const IOobjectList& objects,
39  PtrList<GeometricField<Type, PatchField, GeoMesh> >& fields,
40  const bool syncPar,
41  const bool readOldTime
42 )
43 {
44  typedef GeometricField<Type, PatchField, GeoMesh> GeoField;
45 
46  // Search list of objects for wanted type
47  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
48 
49  const wordList masterNames(fieldNames(fieldObjects, syncPar));
50 
51  fields.setSize(masterNames.size());
52 
53  // Make sure to read in masterNames order.
54 
55  forAll(masterNames, i)
56  {
57  Info<< "Reading " << GeoField::typeName << ' ' << masterNames[i]
58  << endl;
59 
60  const IOobject& io = *fieldObjects[masterNames[i]];
61 
62  fields.set
63  (
64  i,
65  new GeoField
66  (
67  IOobject
68  (
69  io.name(),
70  io.instance(),
71  io.local(),
72  io.db(),
73  IOobject::MUST_READ,
74  IOobject::AUTO_WRITE,
75  io.registerObject()
76  ),
77  mesh,
78  readOldTime
79  )
80  );
81  }
82  return masterNames;
83 }
84 
85 
86 // Read all fields of type. Returns names of fields read. Guarantees all
87 // processors to read fields in same order.
88 template<class GeoField, class Mesh>
90 (
91  const Mesh& mesh,
92  const IOobjectList& objects,
93  PtrList<GeoField>& fields,
94  const bool syncPar
95 )
96 {
97  // Search list of objects for wanted type
98  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
99 
100  const wordList masterNames(fieldNames(fieldObjects, syncPar));
101 
102  fields.setSize(masterNames.size());
103 
104  // Make sure to read in masterNames order.
105 
106  forAll(masterNames, i)
107  {
108  Info<< "Reading " << GeoField::typeName << ' ' << masterNames[i]
109  << endl;
110 
111  const IOobject& io = *fieldObjects[masterNames[i]];
112 
113  fields.set
114  (
115  i,
116  new GeoField
117  (
118  IOobject
119  (
120  io.name(),
121  io.instance(),
122  io.local(),
123  io.db(),
124  IOobject::MUST_READ,
125  IOobject::AUTO_WRITE,
126  io.registerObject()
127  ),
128  mesh
129  )
130  );
131  }
132  return masterNames;
133 }
134 
135 
136 // Read all (non-mesh) fields of type. Returns names of fields read. Guarantees
137 // all processors to read fields in same order.
138 template<class GeoField>
140 (
141  const IOobjectList& objects,
142  PtrList<GeoField>& fields,
143  const bool syncPar
144 )
145 {
146  // Search list of objects for wanted type
147  IOobjectList fieldObjects(objects.lookupClass(GeoField::typeName));
148 
149  const wordList masterNames(fieldNames(fieldObjects, syncPar));
150 
151  fields.setSize(masterNames.size());
152 
153  // Make sure to read in masterNames order.
154 
155  forAll(masterNames, i)
156  {
157  Info<< "Reading " << GeoField::typeName << ' ' << masterNames[i]
158  << endl;
159 
160  const IOobject& io = *fieldObjects[masterNames[i]];
161 
162  fields.set
163  (
164  i,
165  new GeoField
166  (
167  IOobject
168  (
169  io.name(),
170  io.instance(),
171  io.local(),
172  io.db(),
173  IOobject::MUST_READ,
174  IOobject::AUTO_WRITE,
175  io.registerObject()
176  )
177  )
178  );
179  }
180  return masterNames;
181 }
182 
183 
184 template<class GeoField>
185 void Foam::ReadFields
186 (
187  const word& fieldName,
188  const typename GeoField::Mesh& mesh,
189  const wordList& timeNames,
190  objectRegistry& fieldsCache
191 )
192 {
193  // Collect all times that are no longer used
194  {
195  HashSet<word> usedTimes(timeNames);
196 
197  DynamicList<word> unusedTimes(fieldsCache.size());
198 
199  forAllIter(objectRegistry, fieldsCache, timeIter)
200  {
201  const word& tm = timeIter.key();
202  if (!usedTimes.found(tm))
203  {
204  unusedTimes.append(tm);
205  }
206  }
207 
208  //Info<< "Unloading times " << unusedTimes << endl;
209 
210  forAll(unusedTimes, i)
211  {
212  objectRegistry& timeCache = const_cast<objectRegistry&>
213  (
214  fieldsCache.lookupObject<objectRegistry>(unusedTimes[i])
215  );
216  fieldsCache.checkOut(timeCache);
217  }
218  }
219 
220 
221  // Load any new fields
222  forAll(timeNames, i)
223  {
224  const word& tm = timeNames[i];
225 
226  // Create if not found
227  if (!fieldsCache.found(tm))
228  {
229  //Info<< "Creating registry for time " << tm << endl;
230 
231  // Create objectRegistry if not found
232  objectRegistry* timeCachePtr = new objectRegistry
233  (
234  IOobject
235  (
236  tm,
237  tm,
238  fieldsCache,
239  IOobject::NO_READ,
240  IOobject::NO_WRITE
241  )
242  );
243  timeCachePtr->store();
244  }
245 
246  // Obtain cache for current time
247  const objectRegistry& timeCache =
248  fieldsCache.lookupObject<objectRegistry>
249  (
250  tm
251  );
252 
253  // Store field if not found
254  if (!timeCache.found(fieldName))
255  {
256  //Info<< "Loading field " << fieldName
257  // << " for time " << tm << endl;
258 
259  GeoField loadedFld
260  (
261  IOobject
262  (
263  fieldName,
264  tm,
265  mesh.thisDb(),
266  IOobject::MUST_READ,
267  IOobject::NO_WRITE,
268  false
269  ),
270  mesh
271  );
272 
273  // Transfer to timeCache (new objectRegistry and store flag)
274  GeoField* fldPtr = new GeoField
275  (
276  IOobject
277  (
278  fieldName,
279  tm,
280  timeCache,
281  IOobject::NO_READ,
282  IOobject::NO_WRITE
283  ),
284  loadedFld
285  );
286  fldPtr->store();
287  }
288  }
289 }
290 
291 
292 template<class GeoField>
293 void Foam::ReadFields
294 (
295  const word& fieldName,
296  const typename GeoField::Mesh& mesh,
297  const wordList& timeNames,
298  const word& registryName
299 )
300 {
301  ReadFields<GeoField>
302  (
303  fieldName,
304  mesh,
305  timeNames,
306  const_cast<objectRegistry&>
307  (
308  mesh.thisDb().subRegistry(registryName, true)
309  )
310  );
311 }
312 
313 
314 // ************************************************************************* //
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
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
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.
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::Info
messageStream Info
HashSet.H
fieldNames
static List< word > fieldNames
Definition: globalFoam.H:46
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
ReadFields.H
Helper routine to read fields.
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