faMeshReconstructor.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) 2021-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 #include "faMeshReconstructor.H"
29 #include "globalIndex.H"
30 #include "globalMeshData.H"
31 #include "edgeHashes.H"
32 #include "ignoreFaPatch.H"
33 #include "Time.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::faMeshReconstructor::calcAddressing
43 (
44  const labelUList& fvFaceProcAddr
45 )
46 {
47  const globalIndex globalFaceNum(procMesh_.nFaces());
48 
49  // ----------------------
50  // boundaryProcAddressing
51  //
52  // Trivial since non-processor boundary ordering is identical
53 
54  const label nPatches = procMesh_.boundary().size();
55 
56  faBoundaryProcAddr_ = identity(nPatches);
57 
58  // Mark processor patches
59  for
60  (
61  label patchi = procMesh_.boundary().nNonProcessor();
62  patchi < nPatches;
63  ++patchi
64  )
65  {
66  faBoundaryProcAddr_[patchi] = -1;
67  }
68 
69 
70  // ------------------
71  // faceProcAddressing
72  //
73  // Transcribe/rewrite based on finiteVolume faceProcAddressing
74 
75  faFaceProcAddr_ = procMesh_.faceLabels();
76 
77  // From local faceLabels to proc values
78  for (label& facei : faFaceProcAddr_)
79  {
80  // Use finiteVolume info, ignoring face flips
81  facei = mag(fvFaceProcAddr[facei])-1;
82  }
83 
84 
85  // Make global consistent.
86  // Starting at '0', without gaps in numbering etc.
87  // - the sorted order is the obvious solution
88  {
89  globalFaceNum.gather(faFaceProcAddr_, singlePatchFaceLabels_);
90 
91  const labelList globalOrder(Foam::sortedOrder(singlePatchFaceLabels_));
92 
93  singlePatchFaceLabels_ =
94  labelList(singlePatchFaceLabels_, globalOrder);
95 
96  // Set first face to be zero relative to the finiteArea patch
97  // ie, local-face numbering with the first being 0 on any given patch
98  {
99  label patchFirstMeshfacei
100  (
101  singlePatchFaceLabels_.empty()
102  ? 0
103  : singlePatchFaceLabels_.first()
104  );
105  Pstream::broadcast(patchFirstMeshfacei);
106 
107  for (label& facei : faFaceProcAddr_)
108  {
109  facei -= patchFirstMeshfacei;
110  }
111  }
112 
113  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
114 
115  if (Pstream::master())
116  {
117  // Determine the respective local portions of the global ordering
118 
119  labelList procTargets(globalFaceNum.totalSize());
120 
121  for (const label proci : Pstream::allProcs())
122  {
124  (
125  globalFaceNum.range(proci),
126  procTargets
127  ) = proci;
128  }
129 
130  labelList procStarts(globalFaceNum.offsets());
131  labelList procOrders(globalFaceNum.totalSize());
132 
133  for (const label globali : globalOrder)
134  {
135  const label proci = procTargets[globali];
136 
137  procOrders[procStarts[proci]++] =
138  (globali - globalFaceNum.localStart(proci));
139  }
140 
141  // Send the local portions
142  for (const int proci : Pstream::subProcs())
143  {
144  SubList<label> localOrder
145  (
146  procOrders,
147  globalFaceNum.range(proci)
148  );
149 
150  UOPstream toProc(proci, pBufs);
151  toProc << localOrder;
152  }
153 
154  SubList<label> localOrder(procOrders, globalFaceNum.range(0));
155 
156  faFaceProcAddr_ = labelList(faFaceProcAddr_, localOrder);
157  }
158 
159  pBufs.finishedScatters();
160 
161  if (!Pstream::master())
162  {
163  labelList localOrder;
164 
165  UIPstream fromProc(Pstream::masterNo(), pBufs);
166  fromProc >> localOrder;
167 
168  faFaceProcAddr_ = labelList(faFaceProcAddr_, localOrder);
169  }
170  }
171 
172  // Broadcast the same information everywhere
173  Pstream::broadcast(singlePatchFaceLabels_);
174 
175 
176  // ------------------
177  // edgeProcAddressing
178  //
179  // This is more complicated.
180 
181  // Create a single (serial) patch of the finiteArea mesh without a
182  // corresponding volume mesh, otherwise it would be the same as
183  // reconstructPar on the volume mesh (big, slow, etc).
184  //
185  // Do manual point-merge and relabeling to avoid dependency on
186  // finiteVolume pointProcAddressing
187 
188  faPointProcAddr_.clear(); // Final size == procMesh_.nPoints()
189 
190  // 1.
191  // - Topological point merge of the area meshes
192  // - use the local patch faces and points
193 
194  // 2.
195  // - build a single (serial) patch of the finiteArea mesh only
196  // without any point support from the volume mesh
197  // - it may be possible to skip this step, but not obvious how
198 
199  // The collected faces are contiguous per processor, which gives us
200  // the possibility to easily identify their source by using the
201  // global face indexing (if needed).
202  // The ultimate face locations are handled with a separate ordering
203  // list.
204 
205  const uindirectPrimitivePatch& procPatch = procMesh_.patch();
206 
207 
208  {
209  faceList singlePatchProcFaces; // [proc0faces, proc1faces ...]
210  labelList uniqueMeshPointLabels;
211 
212  // Local face points to compact merged point index
213  labelList pointToGlobal;
214 
215  autoPtr<globalIndex> globalPointsPtr =
216  procMesh_.mesh().globalData().mergePoints
217  (
218  procPatch.meshPoints(),
219  procPatch.meshPointMap(), // unused
220  pointToGlobal,
221  uniqueMeshPointLabels
222  );
223 
224  // Gather faces, renumbered for the *merged* points
225  faceList tmpFaces(globalFaceNum.localSize());
226 
227  forAll(tmpFaces, facei)
228  {
229  tmpFaces[facei] =
230  face(pointToGlobal, procPatch.localFaces()[facei]);
231  }
232 
233  globalFaceNum.gather
234  (
235  tmpFaces,
236  singlePatchProcFaces,
239  );
240 
241  globalPointsPtr().gather
242  (
243  UIndirectList<point>
244  (
245  procPatch.points(), // NB: mesh points (non-local)
246  uniqueMeshPointLabels
247  ),
248  singlePatchPoints_
249  );
250 
251  // Transcribe into final assembled order
252  singlePatchFaces_.resize(singlePatchProcFaces.size());
253 
254  // Target face ordering
255  labelList singlePatchProcAddr;
256  globalFaceNum.gather(faFaceProcAddr_, singlePatchProcAddr);
257 
258  forAll(singlePatchProcAddr, facei)
259  {
260  const label targetFacei = singlePatchProcAddr[facei];
261  singlePatchFaces_[targetFacei] = singlePatchProcFaces[facei];
262  }
263 
264  // Use initial equivalent "global" (serial) patch
265  // to establish the correct point and face walking order
266  //
267  // - only meaningful on master
268  const primitivePatch initialPatch
269  (
270  SubList<face>(singlePatchFaces_),
271  singlePatchPoints_
272  );
273 
274  // Grab the faces/points in local point ordering
275  tmpFaces = initialPatch.localFaces();
276  pointField tmpPoints(initialPatch.localPoints());
277 
278  // Equivalent to a flattened meshPointMap
279  labelList mpm(initialPatch.points().size(), -1);
280  {
281  const labelList& mp = initialPatch.meshPoints();
282  forAll(mp, i)
283  {
284  mpm[mp[i]] = i;
285  }
286  }
287  Pstream::broadcast(mpm);
288 
289  // Rewrite pointToGlobal according to the correct point order
290  for (label& pointi : pointToGlobal)
291  {
292  pointi = mpm[pointi];
293  }
294 
295  // Finalize. overwrite
296  faPointProcAddr_ = std::move(pointToGlobal);
297 
298  singlePatchFaces_ = std::move(tmpFaces);
299  singlePatchPoints_ = std::move(tmpPoints);
300  }
301 
302  // Broadcast the same information everywhere
303  Pstream::broadcast(singlePatchFaces_);
304  Pstream::broadcast(singlePatchPoints_);
305 
306  // Now have enough global information to determine global edge mappings
307 
308  // Equivalent "global" (serial) patch
309  const primitivePatch onePatch
310  (
311  SubList<face>(singlePatchFaces_),
312  singlePatchPoints_
313  );
314 
315  faEdgeProcAddr_.clear();
316  faEdgeProcAddr_.resize(procMesh_.nEdges(), -1);
317 
318  {
319  EdgeMap<label> globalEdgeMapping(2*onePatch.nEdges());
320 
321  // Pass 1: edge-hash lookup with edges in "natural" patch order
322 
323  // Can use local edges() directly without remap via meshPoints()
324  // since the previous sorting means that we are already working
325  // with faces that are in the local point order and even
326  // the points themselves are also in the local point order
327  for (const edge& e : onePatch.edges())
328  {
329  globalEdgeMapping.insert(e, globalEdgeMapping.size());
330  }
331 
332  // Lookup proc local edges (in natural order) to global equivalent
333  for (label edgei = 0; edgei < procPatch.nEdges(); ++edgei)
334  {
335  const edge globalEdge(faPointProcAddr_, procPatch.edges()[edgei]);
336 
337  const auto fnd = globalEdgeMapping.cfind(globalEdge);
338 
339  if (fnd.good())
340  {
341  faEdgeProcAddr_[edgei] = fnd.val();
342  }
343  else
344  {
346  << "Failed to find edge " << globalEdge
347  << " this indicates a programming error" << nl
348  << exit(FatalError);
349  }
350  }
351  }
352 
353  // Now have consistent global edge
354  // This of course would all be too easy.
355  // The finiteArea edge lists have been defined as their own sorted
356  // order with sublists etc.
357 
358  // Gather edge ids for nonProcessor boundaries.
359  // These will also be in the serial geometry
360 
361  Map<label> remapGlobal(2*onePatch.nEdges());
362  for (label edgei = 0; edgei < onePatch.nInternalEdges(); ++edgei)
363  {
364  remapGlobal.insert(edgei, remapGlobal.size());
365  }
366 
367  //
368  singlePatchEdgeLabels_.clear();
369  singlePatchEdgeLabels_.resize(procMesh_.boundary().nNonProcessor());
370 
371  forAll(singlePatchEdgeLabels_, patchi)
372  {
373  const faPatch& fap = procMesh_.boundary()[patchi];
374 
375  if (isA<ignoreFaPatch>(fap))
376  {
377  // These are not real edges
378  continue;
379  }
380 
381  labelList& patchEdgeLabels = singlePatchEdgeLabels_[patchi];
382  patchEdgeLabels = fap.edgeLabels();
383 
384  // Renumber from local edges to global edges (natural order)
385  for (label& edgeId : patchEdgeLabels)
386  {
387  edgeId = faEdgeProcAddr_[edgeId];
388  }
389 
390  // OR patchEdgeLabels =
391  // UIndirectList<label>(faEdgeProcAddr_, fap.edgeLabels());
392 
393  // Combine from all processors
394  Pstream::combineReduce(patchEdgeLabels, ListOps::appendEqOp<label>());
395 
396  // Sorted order will be the original non-decomposed order
397  Foam::sort(patchEdgeLabels);
398 
399  for (const label sortedEdgei : patchEdgeLabels)
400  {
401  remapGlobal.insert(sortedEdgei, remapGlobal.size());
402  }
403  }
404 
405  {
406  // Use the map to rewrite the local faEdgeProcAddr_
407 
408  labelList newEdgeProcAddr(faEdgeProcAddr_);
409 
410  label edgei = procMesh_.nInternalEdges();
411 
412  for (const faPatch& fap : procMesh_.boundary())
413  {
414  for (const label patchEdgei : fap.edgeLabels())
415  {
416  const label globalEdgei = faEdgeProcAddr_[patchEdgei];
417 
418  const auto fnd = remapGlobal.cfind(globalEdgei);
419  if (fnd.good())
420  {
421  newEdgeProcAddr[edgei] = fnd.val();
422  ++edgei;
423  }
424  else
425  {
427  << "Failed to find edge " << globalEdgei
428  << " this indicates a programming error" << nl
429  << exit(FatalError);
430  }
431  }
432  }
433 
434  faEdgeProcAddr_ = std::move(newEdgeProcAddr);
435  }
436 }
437 
438 
439 void Foam::faMeshReconstructor::initPatch() const
440 {
441  serialPatchPtr_.reset
442  (
443  new primitivePatch
444  (
445  SubList<face>(singlePatchFaces_),
446  singlePatchPoints_
447  )
448  );
449 }
450 
451 
452 void Foam::faMeshReconstructor::createMesh()
453 {
454  const Time& runTime = procMesh_.thisDb().time();
455 
456  // Time for non-parallel case (w/o functionObjects or libs)
457  serialRunTime_ = Time::NewGlobalTime(runTime);
458 
459 
460  // Trivial polyMesh only containing points and faces.
461  // This is valid, provided we don't use it for anything complicated.
462 
463  serialVolMesh_.reset
464  (
465  new polyMesh
466  (
467  IOobject
468  (
469  procMesh_.mesh().name(), // Volume region name
470  procMesh_.mesh().facesInstance(),
471 
472  *serialRunTime_,
473 
477  ),
478  pointField(singlePatchPoints_), // copy
479  faceList(singlePatchFaces_), // copy
480  labelList(singlePatchFaces_.size(), Zero), // faceOwner
481  labelList(), // faceNeighbour
482  false // no syncPar!
483  )
484  );
485 
486  // A simple area-mesh with one-to-one mapping of faceLabels
487  serialAreaMesh_.reset
488  (
489  new faMesh
490  (
491  *serialVolMesh_,
492  identity(singlePatchFaces_.size()) // faceLabels
493  )
494  );
495 
496  auto& completeMesh = *serialAreaMesh_;
497 
498  // Add in non-processor boundary patches
499  faPatchList completePatches(singlePatchEdgeLabels_.size());
500  label nPatches = 0;
501  forAll(completePatches, patchi)
502  {
503  const labelList& patchEdgeLabels = singlePatchEdgeLabels_[patchi];
504 
505  const faPatch& fap = procMesh_.boundary()[patchi];
506 
507  if (isA<ignoreFaPatch>(fap))
508  {
509  // These are not real edges
510  continue;
511  }
512 
513  const label neiPolyPatchId = fap.ngbPolyPatchIndex();
514 
515  completePatches.set
516  (
517  nPatches,
518  fap.clone
519  (
520  completeMesh.boundary(),
521  patchEdgeLabels,
522  nPatches, // index
523  neiPolyPatchId
524  )
525  );
526 
527  ++nPatches;
528  }
529 
530  completePatches.resize(nPatches);
531 
532  // Serial mesh - no parallel communication
533 
534  const bool oldParRun = Pstream::parRun(false);
535 
536  completeMesh.addFaPatches(completePatches);
537 
538  Pstream::parRun(oldParRun); // Restore parallel state
539 }
540 
541 
542 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
543 
544 Foam::faMeshReconstructor::faMeshReconstructor
545 (
546  const faMesh& procMesh,
547  IOobjectOption::readOption readVolAddressing
548 )
549 :
550  procMesh_(procMesh),
551  errors_(0)
552 {
553  if (!Pstream::parRun())
554  {
556  << "Can only be called in parallel!!" << nl
557  << exit(FatalError);
558  }
559 
560  IOobject ioAddr
561  (
562  "faceProcAddressing",
563  procMesh_.mesh().facesInstance(),
564 
565  // Or search?
566  // procMesh_.time().findInstance
567  // (
568  // // Search for polyMesh face instance
569  // // mesh.facesInstance()
570  // procMesh_.mesh().meshDir(),
571  // "faceProcAddressing"
572  // ),
573 
575  procMesh_.mesh(), // The polyMesh db
576 
577  readVolAddressing, // Read option
580  );
581 
582  // Require faceProcAddressing from finiteVolume decomposition
583  labelIOList fvFaceProcAddr(ioAddr);
584 
585  // Check if any/all where read.
586  // Use 'headerClassName' for checking
587  bool fileOk
588  (
589  fvFaceProcAddr.isAnyRead()
590  && fvFaceProcAddr.isHeaderClass<labelIOList>()
591  );
592 
593  Pstream::reduceAnd(fileOk);
594 
595  if (fileOk)
596  {
597  calcAddressing(fvFaceProcAddr);
598  }
599  else
600  {
601  errors_ = 1;
602  }
603 }
604 
605 
606 Foam::faMeshReconstructor::faMeshReconstructor
607 (
608  const faMesh& procMesh,
609  const labelUList& fvFaceProcAddressing
610 )
611 :
612  procMesh_(procMesh),
613  errors_(0)
614 {
615  if (!Pstream::parRun())
616  {
618  << "Can only be called in parallel!!" << nl
619  << exit(FatalError);
620  }
622  calcAddressing(fvFaceProcAddressing);
623 }
624 
625 
626 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
629 {}
630 
631 
632 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
633 
635 {
636  serialAreaMesh_.reset(nullptr);
637  serialVolMesh_.reset(nullptr);
638  serialRunTime_.reset(nullptr);
639 }
640 
641 
643 {
644  if (!serialPatchPtr_)
645  {
646  initPatch();
647  }
648 
649  return *serialPatchPtr_;
650 }
651 
652 
654 {
655  if (!serialPatchPtr_)
656  {
657  initPatch();
658  }
659 
660  return *serialPatchPtr_;
661 }
662 
663 
665 {
666  if (!serialAreaMesh_)
667  {
668  const_cast<faMeshReconstructor&>(*this).createMesh();
669  }
670 
671  return *serialAreaMesh_;
672 }
673 
676 {
677  writeAddressing(procMesh_.mesh().facesInstance());
678 }
679 
680 
682 {
683  // Write copies
684 
685  IOobject ioAddr
686  (
687  "procAddressing",
688  timeName,
690  procMesh_.thisDb(),
694  );
695 
696  // boundaryProcAddressing
697  ioAddr.rename("boundaryProcAddressing");
698  IOListRef<label>(ioAddr, faBoundaryProcAddr_).write();
699 
700  // faceProcAddressing
701  ioAddr.rename("faceProcAddressing");
702  IOListRef<label>(ioAddr, faFaceProcAddr_).write();
703 
704  // pointProcAddressing
705  ioAddr.rename("pointProcAddressing");
706  IOListRef<label>(ioAddr, faPointProcAddr_).write();
708  // edgeProcAddressing
709  ioAddr.rename("edgeProcAddressing");
710  IOListRef<label>(ioAddr, faEdgeProcAddr_).write();
711 }
712 
715 {
716  writeMesh(procMesh_.mesh().facesInstance());
717 }
718 
719 
721 {
722  const faMesh& fullMesh = this->mesh();
723 
725 
726  auto oldHandler = fileOperation::fileHandler(writeHandler);
727  const bool oldDistributed = fileHandler().distributed(true);
728 
729  if (UPstream::master())
730  {
731  const bool oldParRun = UPstream::parRun(false);
732 
733  IOobject io(fullMesh.boundary());
734 
735  io.rename("faceLabels");
736  IOListRef<label>(io, singlePatchFaceLabels_).write();
737 
738  fullMesh.boundary().write();
739 
740  UPstream::parRun(oldParRun);
741  }
742 
743  // Restore settings
744  fileHandler().distributed(oldDistributed);
745  (void) fileOperation::fileHandler(oldHandler);
746 }
747 
748 
749 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:59
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
void writeMesh() const
Write equivalent mesh information at the polyMesh faceInstances time.
const faMesh & mesh() const
Serial equivalent faMesh.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static autoPtr< Time > NewGlobalTime()
Construct (dummy) global Time - no functionObjects or libraries, using the global path information st...
Definition: TimeNew.C:78
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
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
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:1176
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
virtual void rename(const word &newName)
Rename the object.
Definition: IOobject.H:677
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
engineTime & runTime
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
static int debug
Debug flag.
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
const faGlobalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: faMesh.C:1033
SubList< label > subList
Declare type of subList.
Definition: List.H:144
Ignore writing from objectRegistry::writeObject()
A bare-bones reconstructor for finiteArea meshes when processor meshes are available (in parallel) bu...
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
static const faMesh & mesh(const polyMesh &pMesh)
The single-region finite-area region on the polyMesh. Uses lookupObject semantics - Fatal if non-exis...
Definition: faMesh.C:74
const primitivePatch & patch() const
Serial equivalent patch.
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:61
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
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...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:109
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
word timeName
Definition: getTimeIndex.H:3
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:103
static void reduceAnd(bool &value, const label communicator=worldComm)
Logical (and) reduction (MPI_AllReduce)
"scheduled" : (MPI_Send, MPI_Recv)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:484
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:31
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
label nEdges() const noexcept
Number of local mesh edges.
Definition: faMeshI.H:43
void writeAddressing() const
Write proc addressing at the polyMesh faceInstances time.
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:701
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:803
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:49
label nNonProcessor() const
The number of patches before the first processor patch.
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.
"nonBlocking" : (MPI_Isend, MPI_Irecv)
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:1185
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A IOList wrapper for writing external data.
Definition: IOList.H:151
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const objectRegistry & thisDb() const noexcept
Return the object registry.
Do not request registration (bool: false)
const dimensionedScalar mp
Proton mass.
static autoPtr< fileOperation > NewUncollated()
The commonly used uncollatedFileOperation.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
readOption
Enumeration defining read preferences.