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-2023 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  const UPtrList<const GeoField> other
171  (
172  mesh.thisDb().objectRegistry::template cobjects<GeoField>()
173  );
174 
175  for (const GeoField& field : other)
176  {
177  if (!field.ownedByRegistry())
178  {
179  const_cast<GeoField&>(field).checkOut();
180  }
181  }
182  }
183 
184  // Early exit
185  return;
186  }
187 
188 
189  // We are decomposing if none of the subprocs has a mesh
190  bool decompose = true;
191  if (UPstream::master())
192  {
193  for (const int proci : UPstream::subProcs())
194  {
195  if (haveMeshOnProc.test(proci))
196  {
197  decompose = false;
198  break;
199  }
200  }
201  }
202  Pstream::broadcast(decompose);
203 
204 
205  if (decompose && UPstream::master())
206  {
207  const bool oldParRun = UPstream::parRun(false);
208 
209  forAll(masterNames, i)
210  {
211  const word& name = masterNames[i];
212  IOobject& io = *objects[name];
214 
215  // Load field (but not oldTime)
216  readField(io, mesh, i, fields);
217  }
218 
219  UPstream::parRun(oldParRun);
220  }
221  else if
222  (
223  !decompose
224  &&
225  // Has read-handler : use it to decide if reading is possible
226  // No read-handler : decide based on the presence of a mesh
227  (
228  readHandlerPtr
229  ? readHandlerPtr->good()
230  : haveMeshOnProc.test(UPstream::myProcNo())
231  )
232  )
233  {
234  const label oldWorldComm = UPstream::worldComm;
235  refPtr<fileOperation> oldHandler;
236 
237  if (readHandlerPtr)
238  {
239  // Swap read fileHandler for read fields
240  oldHandler = fileOperation::fileHandler(*readHandlerPtr);
242  }
243 
244  forAll(masterNames, i)
245  {
246  const word& name = masterNames[i];
247  IOobject& io = *objects[name];
249 
250  // Load field (but not oldTime)
251  readField(io, mesh, i, fields);
252  }
253 
254  if (readHandlerPtr)
255  {
256  // Restore fileHandler
257  *readHandlerPtr = fileOperation::fileHandler(oldHandler);
258  UPstream::commWorld(oldWorldComm);
259  }
260  }
261 
262 
263  // Missing fields on any processors?
264  // - construct from dictionary
265 
266  PtrList<dictionary> fieldDicts;
267 
268  if (UPstream::master())
269  {
270  // Broadcast zero sized fields everywhere (if needed)
271  // Send like a list of dictionaries
272 
273  OPBstream toProcs(UPstream::masterNo()); // worldComm
274 
275  const label nDicts = (subsetter ? fields.size() : label(0));
276 
277  toProcs << nDicts << token::BEGIN_LIST; // Begin list
278 
279  if (nDicts && subsetter)
280  {
281  // Disable communication for interpolate() method
282  const bool oldParRun = UPstream::parRun(false);
283 
284  for (const auto& fld : fields)
285  {
286  tmp<GeoField> tsubfld = subsetter->interpolate(fld);
287 
288  // Surround each with {} as dictionary entry
289  toProcs.beginBlock();
290  toProcs << tsubfld();
291  toProcs.endBlock();
292  }
293 
294  UPstream::parRun(oldParRun); // Restore state
295  }
296 
297  toProcs << token::END_LIST << token::NL; // End list
298  }
299  else
300  {
301  // Receive the broadcast...
302  IPBstream fromMaster(UPstream::masterNo()); // worldComm
303 
304  // But only consume where needed...
305  if (!haveMeshOnProc.test(UPstream::myProcNo()))
306  {
307  fromMaster >> fieldDicts;
308  }
309  }
310 
311 
312  // Use the received dictionaries (if any) to create missing fields.
313 
314  // Disable communication when constructing from dictionary
315  const bool oldParRun = UPstream::parRun(false);
316 
317  IOobject noreadIO
318  (
319  "none",
320  mesh.time().timeName(),
321  mesh.thisDb(),
325  );
326 
327  forAll(fieldDicts, i)
328  {
329  noreadIO.resetHeader(masterNames[i]);
330 
331  fields.emplace_set(i, noreadIO, mesh, fieldDicts[i]);
332  }
333 
334  UPstream::parRun(oldParRun); // Restore any changes
335 
336 
337  // Finally. Can checkOut of registry as required
338  if (deregister)
339  {
341  for (auto& fld : fields)
342  {
344 
345  // Ensure it is not destroyed by polyMesh deletion
346  fld.checkOut();
347  }
349 
350  // Extra safety - remove all such types
351  HashTable<const GeoField*> other
352  (
353  mesh.thisDb().objectRegistry::template lookupClass<GeoField>()
354  );
355 
356  forAllConstIters(other, iter)
357  {
358  GeoField& fld = const_cast<GeoField&>(*iter.val());
359 
360  if (!fld.ownedByRegistry())
361  {
362  fld.checkOut();
363  }
364  }
365  }
366 }
367 
368 
369 template<class GeoField, class MeshSubsetter>
371 (
372  const bitSet& haveMeshOnProc,
373  const MeshSubsetter* subsetter,
374  const typename GeoField::Mesh& mesh,
375  IOobjectList& allObjects,
376  PtrList<GeoField>& fields,
377  const bool deregister
378 )
379 {
380  readFieldsImpl
381  (
382  nullptr, // readHandler
383  haveMeshOnProc,
384  subsetter,
385 
386  mesh,
387  allObjects,
388  fields,
389  deregister
390  );
391 }
392 
393 
394 template<class GeoField, class MeshSubsetter>
396 (
397  const boolUList& haveMeshOnProc,
398  const MeshSubsetter* subsetter,
399  const typename GeoField::Mesh& mesh,
400  IOobjectList& allObjects,
402  const bool deregister
403 )
404 {
405  readFieldsImpl
406  (
407  nullptr, // readHandler
408  haveMeshOnProc,
409  subsetter,
410 
411  mesh,
412  allObjects,
413  fields,
414  deregister
415  );
416 }
417 
418 
419 template<class GeoField, class MeshSubsetter>
421 (
422  const boolUList& haveMeshOnProc,
423  const typename GeoField::Mesh& mesh,
424  const autoPtr<MeshSubsetter>& subsetter,
425  IOobjectList& allObjects,
427  const bool deregister
428 )
429 {
430  readFieldsImpl
431  (
432  nullptr, // readHandler
433  haveMeshOnProc,
434  subsetter.get(),
435 
436  mesh,
437  allObjects,
438  fields,
439  deregister
440  );
441 }
442 
443 
444 template<class GeoField, class MeshSubsetter>
446 (
447  const boolUList& haveMeshOnProc,
448  refPtr<fileOperation>& readHandler,
449  const typename GeoField::Mesh& mesh,
450  const autoPtr<MeshSubsetter>& subsetter,
451  IOobjectList& allObjects,
453  const bool deregister
454 )
455 {
456  readFieldsImpl
457  (
458  &readHandler,
459  haveMeshOnProc,
460  subsetter.get(),
461 
462  mesh,
463  allObjects,
464  fields,
465  deregister
466  );
467 }
468 
469 
470 // ************************************************************************* //
rDeltaTY field()
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:598
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:1049
Begin list [isseparator].
Definition: token.H:161
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
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:1074
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:409
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.
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
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:429
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:1082
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:1185
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:172
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
Request registration (bool: true)
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