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-2024 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 idx,
37 )
38 {
39  fields.emplace_set(idx, io, mesh);
40 }
41 
42 
43 template<class Type, template<class> class PatchField, class GeoMesh>
45 (
46  const IOobject& io,
47  const typename GeoMesh::Mesh& mesh,
48  const label idx,
50 )
51 {
52  fields.emplace_set(idx, io, mesh, false); // readOldTime = false
53 }
54 
55 
56 template<class Type, template<class> class PatchField, class GeoMesh>
58 (
59  const typename GeoMesh::Mesh& mesh,
60  const IOobjectList& objects,
62  const bool readOldTime
63 )
64 {
66 
67  // GeoField fields - sorted for consistent order on all processors
68  const UPtrList<const IOobject> fieldObjects(objects.csorted<GeoField>());
69 
70  // Construct the fields
71  fields.resize_null(fieldObjects.size());
72 
73  forAll(fieldObjects, i)
74  {
75  fields.emplace_set(i, fieldObjects[i], mesh, readOldTime);
76  }
77 }
78 
79 
80 template<class Mesh, class GeoField>
82 (
83  const Mesh& mesh,
84  const IOobjectList& objects,
85  PtrList<GeoField>& fields
86 )
87 {
88  // GeoField fields - sorted for consistent order on all processors
89  const UPtrList<const IOobject> fieldObjects(objects.csorted<GeoField>());
90 
91  // Construct the fields
92  fields.resize_null(fieldObjects.size());
93 
94  forAll(fieldObjects, i)
95  {
96  fields.emplace_set(i, fieldObjects[i], mesh);
97  }
98 }
99 
100 
101 template<class BoolListType, class GeoField, class MeshSubsetter>
102 void Foam::fieldsDistributor::readFieldsImpl
103 (
104  refPtr<fileOperation>* readHandlerPtr, // Can be nullptr
105  const BoolListType& haveMeshOnProc,
106  const MeshSubsetter* subsetter,
107  const typename GeoField::Mesh& mesh,
108  IOobjectList& allObjects,
109  PtrList<GeoField>& fields,
110  const bool deregister
111 )
112 {
113  // Get my objects of type
114  IOobjectList objects(allObjects.lookupClass<GeoField>());
115 
116  // Check that we all have all objects
117  wordList objectNames = objects.sortedNames();
118 
119  // Get master names
120  wordList masterNames(objectNames);
121  Pstream::broadcast(masterNames);
122 
123  if
124  (
125  haveMeshOnProc.test(UPstream::myProcNo())
126  && objectNames != masterNames
127  )
128  {
130  << "Objects not synchronised across processors." << nl
131  << "Master has " << flatOutput(masterNames) << nl
132  << "Processor " << UPstream::myProcNo()
133  << " has " << flatOutput(objectNames)
134  << exit(FatalError);
135  }
136 
137  #ifdef FULLDEBUG
138  {
139  // A bit more checking - this may not be strictly correct.
140  // The master must know about everyone, but the others really
141  // only need to know about themselves.
142  // Can broadcast decompose = yes/no from master
143 
144  bitSet localValues(haveMeshOnProc);
145  bitSet masterValues(localValues);
146  Pstream::broadcast(masterValues);
147 
148  localValues ^= masterValues;
149 
150  if (localValues.any())
151  {
153  << "haveMeshOnProc not synchronised across processors." << nl
154  << "Processor " << UPstream::myProcNo()
155  << " differs at these positions: "
156  << flatOutput(localValues.sortedToc()) << nl
157  << exit(FatalError);
158  }
159  }
160  #endif
161 
162 
163  fields.resize_null(masterNames.size());
164 
165  if (fields.empty())
166  {
167  if (deregister)
168  {
169  // Extra safety - remove all such types
170  for
171  (
172  const GeoField& fld
173  : mesh.thisDb().objectRegistry::template csorted<GeoField>()
174  )
175  {
176  if (!fld.ownedByRegistry())
177  {
178  const_cast<GeoField&>(fld).checkOut();
179  }
180  }
181  }
182 
183  // Early exit
184  return;
185  }
186 
187 
188  // We are decomposing if none of the subprocs has a mesh
189  bool decompose = true;
190  if (UPstream::master())
191  {
192  for (const int proci : UPstream::subProcs())
193  {
194  if (haveMeshOnProc.test(proci))
195  {
196  decompose = false;
197  break;
198  }
199  }
200  }
201  Pstream::broadcast(decompose);
202 
203 
204  if (decompose && UPstream::master())
205  {
206  const bool oldParRun = UPstream::parRun(false);
207 
208  forAll(masterNames, i)
209  {
210  const word& name = masterNames[i];
211  IOobject& io = *objects[name];
213 
214  // Load field (but not oldTime)
215  readField(io, mesh, i, fields);
216  }
217 
218  UPstream::parRun(oldParRun);
219  }
220  else if
221  (
222  !decompose
223  &&
224  // Has read-handler : use it to decide if reading is possible
225  // No read-handler : decide based on the presence of a mesh
226  (
227  readHandlerPtr
228  ? readHandlerPtr->good()
229  : haveMeshOnProc.test(UPstream::myProcNo())
230  )
231  )
232  {
233  const label oldWorldComm = UPstream::worldComm;
234  refPtr<fileOperation> oldHandler;
235 
236  if (readHandlerPtr)
237  {
238  // Swap read fileHandler for read fields
239  oldHandler = fileOperation::fileHandler(*readHandlerPtr);
241  }
242 
243  forAll(masterNames, i)
244  {
245  const word& name = masterNames[i];
246  IOobject& io = *objects[name];
248 
249  // Load field (but not oldTime)
250  readField(io, mesh, i, fields);
251  }
252 
253  if (readHandlerPtr)
254  {
255  // Restore fileHandler
256  *readHandlerPtr = fileOperation::fileHandler(oldHandler);
257  UPstream::commWorld(oldWorldComm);
258  }
259  }
260 
261 
262  // Missing fields on any processors?
263  // - construct from dictionary
264 
265  PtrList<dictionary> fieldDicts;
266 
267  if (UPstream::master())
268  {
269  // Broadcast zero sized fields everywhere (if needed)
270  // Send like a list of dictionaries
271 
272  OPBstream toProcs(UPstream::worldComm);
273 
274  const label nDicts = (subsetter ? fields.size() : label(0));
275 
276  toProcs << nDicts << token::BEGIN_LIST; // Begin list
277 
278  if (nDicts && subsetter)
279  {
280  // Disable communication for interpolate() method
281  const bool oldParRun = UPstream::parRun(false);
282 
283  for (const auto& fld : fields)
284  {
285  tmp<GeoField> tsubfld = subsetter->interpolate(fld);
286 
287  // Surround each with {} as dictionary entry
288  toProcs.beginBlock();
289  toProcs << tsubfld();
290  toProcs.endBlock();
291  }
292 
293  UPstream::parRun(oldParRun); // Restore state
294  }
295 
296  toProcs << token::END_LIST << token::NL; // End list
297  }
298  else
299  {
300  // Receive the broadcast...
301  IPBstream fromMaster(UPstream::worldComm);
302 
303  // But only consume where needed...
304  if (!haveMeshOnProc.test(UPstream::myProcNo()))
305  {
306  fromMaster >> fieldDicts;
307  }
308  }
309 
310 
311  // Use the received dictionaries (if any) to create missing fields.
312 
313  // Disable communication when constructing from dictionary
314  const bool oldParRun = UPstream::parRun(false);
315 
316  IOobject noreadIO
317  (
318  "none",
319  mesh.time().timeName(),
320  mesh.thisDb(),
324  );
325 
326  forAll(fieldDicts, i)
327  {
328  noreadIO.resetHeader(masterNames[i]);
329 
330  fields.emplace_set(i, noreadIO, mesh, fieldDicts[i]);
331  }
332 
333  UPstream::parRun(oldParRun); // Restore any changes
334 
335 
336  // Finally. Can checkOut of registry as required
337  if (deregister)
338  {
340  for (auto& fld : fields)
341  {
343 
344  // Ensure it is not destroyed by polyMesh deletion
345  fld.checkOut();
346  }
348 
349  // Extra safety - remove all such types
350  for
351  (
352  const GeoField& fld
353  : mesh.thisDb().objectRegistry::template csorted<GeoField>()
354  )
355  {
356  if (!fld.ownedByRegistry())
357  {
358  const_cast<GeoField&>(fld).checkOut();
359  }
360  }
361  }
362 }
363 
364 
365 template<class GeoField, class MeshSubsetter>
367 (
368  const bitSet& haveMeshOnProc,
369  const MeshSubsetter* subsetter,
370  const typename GeoField::Mesh& mesh,
371  IOobjectList& allObjects,
372  PtrList<GeoField>& fields,
373  const bool deregister
374 )
375 {
376  readFieldsImpl
377  (
378  nullptr, // readHandler
379  haveMeshOnProc,
380  subsetter,
381 
382  mesh,
383  allObjects,
384  fields,
385  deregister
386  );
387 }
388 
389 
390 template<class GeoField, class MeshSubsetter>
392 (
393  const boolUList& haveMeshOnProc,
394  const MeshSubsetter* subsetter,
395  const typename GeoField::Mesh& mesh,
396  IOobjectList& allObjects,
398  const bool deregister
399 )
400 {
401  readFieldsImpl
402  (
403  nullptr, // readHandler
404  haveMeshOnProc,
405  subsetter,
406 
407  mesh,
408  allObjects,
409  fields,
410  deregister
411  );
412 }
413 
414 
415 template<class GeoField, class MeshSubsetter>
417 (
418  const boolUList& haveMeshOnProc,
419  const typename GeoField::Mesh& mesh,
420  const autoPtr<MeshSubsetter>& subsetter,
421  IOobjectList& allObjects,
423  const bool deregister
424 )
425 {
426  readFieldsImpl
427  (
428  nullptr, // readHandler
429  haveMeshOnProc,
430  subsetter.get(),
431 
432  mesh,
433  allObjects,
434  fields,
435  deregister
436  );
437 }
438 
439 
440 template<class GeoField, class MeshSubsetter>
442 (
443  const boolUList& haveMeshOnProc,
444  refPtr<fileOperation>& readHandler,
445  const typename GeoField::Mesh& mesh,
446  const autoPtr<MeshSubsetter>& subsetter,
447  IOobjectList& allObjects,
449  const bool deregister
450 )
451 {
452  readFieldsImpl
453  (
454  &readHandler,
455  haveMeshOnProc,
456  subsetter.get(),
457 
458  mesh,
459  allObjects,
460  fields,
461  deregister
462  );
463 }
464 
465 
466 // ************************************************************************* //
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:608
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
wordList sortedNames() const
The sorted names of the IOobjects.
Definition: IOobjectList.C:250
Newline [isspace].
Definition: token.H:130
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
Begin list [isseparator].
Definition: token.H:161
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
Generic GeometricField class.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:421
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
writeOption writeOpt() const noexcept
Get the write option.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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.
End list [isseparator].
Definition: token.H:162
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
MESH Mesh
Definition: GeoMesh.H:58
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds)
Definition: UPstream.H:441
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
List of word.
Definition: fileName.H:59
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)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
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
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
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:1197
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
Request registration (bool: true)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225