fieldsDistributorTemplates.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) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 template<class GeoField>
32 (
33  const IOobject& io,
34  const typename GeoField::Mesh& mesh,
35  const label i,
37 )
38 {
39  fields.set(i, new GeoField(io, mesh));
40 }
41 
42 template<class Type, template<class> class PatchField, class GeoMesh>
44 (
45  const IOobject& io,
46  const typename GeoMesh::Mesh& mesh,
47  const label i,
49 )
50 {
51  fields.set
52  (
53  i,
55  );
56 }
57 
58 
59 template<class Type, template<class> class PatchField, class GeoMesh>
61 (
62  const typename GeoMesh::Mesh& mesh,
63  const IOobjectList& objects,
65  const bool readOldTime
66 )
67 {
69 
70  // GeoField fields - sorted for consistent order on all processors
71  UPtrList<const IOobject> fieldObjects(objects.sorted<GeoField>());
72 
73  // Construct the fields
74  fields.resize(fieldObjects.size());
75 
76  forAll(fieldObjects, i)
77  {
78  fields.set(i, new GeoField(fieldObjects[i], mesh, readOldTime));
79  }
80 }
81 
82 
83 template<class Mesh, class GeoField>
85 (
86  const Mesh& mesh,
87  const IOobjectList& objects,
88  PtrList<GeoField>& fields
89 )
90 {
91  // GeoField fields - sorted for consistent order on all processors
92  UPtrList<const IOobject> fieldObjects(objects.sorted<GeoField>());
93 
94  // Construct the fields
95  fields.resize(fieldObjects.size());
96 
97  forAll(fieldObjects, i)
98  {
99  fields.set(i, new GeoField(fieldObjects[i], mesh));
100  }
101 }
102 
103 
104 template<class BoolListType, class GeoField, class MeshSubsetter>
105 void Foam::fieldsDistributor::readFieldsImpl
106 (
107  const BoolListType& haveMeshOnProc,
108  const MeshSubsetter* subsetter,
109  const typename GeoField::Mesh& mesh,
110  IOobjectList& allObjects,
111  PtrList<GeoField>& fields,
112  const bool deregister
113 )
114 {
115  // Get my objects of type
116  IOobjectList objects(allObjects.lookupClass<GeoField>());
117 
118  // Check that we all have all objects
119  wordList objectNames = objects.sortedNames();
120 
121  // Get master names
122  wordList masterNames(objectNames);
123  Pstream::broadcast(masterNames);
124 
125  if (haveMeshOnProc.test(Pstream::myProcNo()) && objectNames != masterNames)
126  {
128  << "Objects not synchronised across processors." << nl
129  << "Master has " << flatOutput(masterNames) << nl
130  << "Processor " << Pstream::myProcNo()
131  << " has " << flatOutput(objectNames)
132  << exit(FatalError);
133  }
134 
135  fields.clear();
136  fields.resize(masterNames.size());
137 
138  if (fields.empty())
139  {
140  if (deregister)
141  {
142  // Extra safety - remove all such types
143  HashTable<const GeoField*> other
144  (
145  mesh.thisDb().objectRegistry::template lookupClass<GeoField>()
146  );
147 
148  forAllConstIters(other, iter)
149  {
150  GeoField& fld = const_cast<GeoField&>(*iter.val());
151 
152  if (!fld.ownedByRegistry())
153  {
154  fld.checkOut();
155  }
156  }
157  }
158 
159  // Early exit
160  return;
161  }
162 
163 
164  // Have master send all fields to processors that don't have a mesh. The
165  // issue is if a patchField does any parallel operations inside its
166  // construct-from-dictionary. This will not work when going to more
167  // processors (e.g. decompose = 1 -> many) ! We could make a special
168  // exception for decomposePar but nicer would be to have read-communicator
169  // ... For now detect if decomposing & disable parRun
170  if (Pstream::master())
171  {
172  // Work out if we're decomposing - none of the subprocs has a mesh
173  bool decompose = true;
174  for (const int proci : Pstream::subProcs())
175  {
176  if (haveMeshOnProc.test(proci))
177  {
178  decompose = false;
179  break;
180  }
181  }
182 
183  const bool oldParRun = Pstream::parRun();
184  if (decompose)
185  {
186  Pstream::parRun(false);
187  }
188 
189  forAll(masterNames, i)
190  {
191  const word& name = masterNames[i];
192  IOobject& io = *objects[name];
194 
195  // Load field (but not oldTime)
196  readField(io, mesh, i, fields);
197  }
198 
199  Pstream::parRun(oldParRun); // Restore any changes
200  }
201  else if (haveMeshOnProc.test(Pstream::myProcNo()))
202  {
203  // Have mesh so just try to load
204  forAll(masterNames, i)
205  {
206  const word& name = masterNames[i];
207  IOobject& io = *objects[name];
209 
211 
212  // Load field (but not oldTime)
213  readField(io, mesh, i, fields);
214  }
215  }
216 
217 
218  // Missing fields on any processors?
219  // - construct from dictionary
220 
221  PtrList<dictionary> fieldDicts;
222 
223  if (Pstream::master())
224  {
225  // Broadcast zero sized fields everywhere (if needed)
226  // Send like a list of dictionaries
227 
228  OPBstream toProcs(UPstream::masterNo()); // worldComm
229 
230  const label nDicts = (subsetter ? fields.size() : label(0));
231 
232  toProcs << nDicts << token::BEGIN_LIST; // Begin list
233 
234  if (nDicts && subsetter)
235  {
236  // Disable communication for interpolate() method
237  const bool oldParRun = Pstream::parRun(false);
238 
239  for (const auto& fld : fields)
240  {
241  tmp<GeoField> tsubfld = subsetter->interpolate(fld);
242 
243  // Surround each with {} as dictionary entry
244  toProcs.beginBlock();
245  toProcs << tsubfld();
246  toProcs.endBlock();
247  }
248 
249  Pstream::parRun(oldParRun); // Restore state
250  }
251 
252  toProcs << token::END_LIST << token::NL; // End list
253  }
254  else
255  {
256  // Receive the broadcast...
257  IPBstream fromMaster(UPstream::masterNo()); // worldComm
258 
259  // But only consume where needed...
260  if (!haveMeshOnProc.test(Pstream::myProcNo()))
261  {
262  fromMaster >> fieldDicts;
263  }
264  }
265 
266 
267  // Use the received dictionaries (if any) to create missing fields.
268 
269  // Disable communication when constructing from dictionary
270  const bool oldParRun = Pstream::parRun(false);
271 
272  forAll(fieldDicts, i)
273  {
274  fields.set
275  (
276  i,
277  new GeoField
278  (
279  IOobject
280  (
281  masterNames[i],
282  mesh.time().timeName(),
283  mesh.thisDb(),
286  ),
287  mesh,
288  fieldDicts[i]
289  )
290  );
291  }
292 
293  Pstream::parRun(oldParRun); // Restore any changes
294 
295 
296  // Finally. Can checkOut of registry as required
297  if (deregister)
298  {
300  for (auto& fld : fields)
301  {
303 
304  // Ensure it is not destroyed by polyMesh deletion
305  fld.checkOut();
306  }
308 
309  // Extra safety - remove all such types
310  HashTable<const GeoField*> other
311  (
312  mesh.thisDb().objectRegistry::template lookupClass<GeoField>()
313  );
314 
315  forAllConstIters(other, iter)
316  {
317  GeoField& fld = const_cast<GeoField&>(*iter.val());
318 
319  if (!fld.ownedByRegistry())
320  {
321  fld.checkOut();
322  }
323  }
324  }
325 }
326 
327 
328 template<class GeoField, class MeshSubsetter>
330 (
331  const bitSet& haveMeshOnProc,
332  const MeshSubsetter* subsetter,
333  const typename GeoField::Mesh& mesh,
334  IOobjectList& allObjects,
335  PtrList<GeoField>& fields,
336  const bool deregister
337 )
338 {
339  readFieldsImpl
340  (
341  haveMeshOnProc,
342  subsetter,
343 
344  mesh,
345  allObjects,
346  fields,
347  deregister
348  );
349 }
350 
351 
352 template<class GeoField, class MeshSubsetter>
354 (
355  const boolList& haveMeshOnProc,
356  const MeshSubsetter* subsetter,
357  const typename GeoField::Mesh& mesh,
358  IOobjectList& allObjects,
360  const bool deregister
361 )
362 {
363  readFieldsImpl
364  (
365  haveMeshOnProc,
366  subsetter,
367 
368  mesh,
369  allObjects,
370  fields,
371  deregister
372  );
373 }
374 
375 
376 template<class GeoField, class MeshSubsetter>
378 (
379  const boolList& haveMeshOnProc,
380  const typename GeoField::Mesh& mesh,
381  const autoPtr<MeshSubsetter>& subsetter,
382  IOobjectList& allObjects,
384  const bool deregister
385 )
386 {
387  readFieldsImpl
388  (
389  haveMeshOnProc,
390  subsetter.get(),
391 
392  mesh,
393  allObjects,
394  fields,
395  deregister
396  );
397 }
398 
399 
400 // ************************************************************************* //
static void readField(const IOobject &io, const typename GeoField::Mesh &mesh, const label i, PtrList< GeoField > &fields)
Generic mesh-based field reading.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
UPtrList< const IOobject > sorted() const
The sorted list of IOobjects.
Definition: IOobjectList.C:254
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:301
Newline [isspace].
Definition: token.H:127
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:639
Begin list [isseparator].
Definition: token.H:158
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all processes in communicator.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:377
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
writeOption writeOpt() const noexcept
Get the write option.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
static void readFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool readOldTime)
Read fields and store on the pointer list.
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:664
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
End list [isseparator].
Definition: token.H:159
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:100
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
MESH Mesh
Definition: GeoMesh.H:58
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
List< word > wordList
A List of words.
Definition: fileName.H:58
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static bool master(const label communicator=worldComm)
Am I the master rank.
Definition: UPstream.H:672
Nothing to be read.
Automatically write from objectRegistry::writeObject()
T * get() noexcept
Return pointer to managed object without nullptr checking.
Definition: autoPtr.H:216
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:757
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225