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 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 "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::scatter(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.size());
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.size());
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.finishedSends();
160 
161  if (!Pstream::master())
162  {
163  labelList localOrder;
164 
165  UIPstream fromProc(Pstream::master(), pBufs);
166  fromProc >> localOrder;
167 
168  faFaceProcAddr_ = labelList(faFaceProcAddr_, localOrder);
169  }
170  }
171 
172  // Broadcast the same information everywhere
173  Pstream::scatter(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::scatter(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::scatter(singlePatchFaces_);
304  Pstream::scatter(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.found())
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  labelList& patchEdgeLabels = singlePatchEdgeLabels_[patchi];
376  patchEdgeLabels = fap.edgeLabels();
377 
378  // Renumber from local edges to global edges (natural order)
379  for (label& edgeId : patchEdgeLabels)
380  {
381  edgeId = faEdgeProcAddr_[edgeId];
382  }
383 
384  // OR patchEdgeLabels =
385  // UIndirectList<label>(faEdgeProcAddr_, fap.edgeLabels());
386 
387 
388  // Collect from all processors
390  (
391  patchEdgeLabels,
392  ListOps::appendEqOp<label>()
393  );
394 
395  // Sorted order will be the original non-decomposed order
396  Foam::sort(patchEdgeLabels);
397 
398  for (const label sortedEdgei : patchEdgeLabels)
399  {
400  remapGlobal.insert(sortedEdgei, remapGlobal.size());
401  }
402  }
403 
404  {
405  // Use the map to rewrite the local faEdgeProcAddr_
406 
407  labelList newEdgeProcAddr(faEdgeProcAddr_);
408 
409  label edgei = procMesh_.nInternalEdges();
410 
411  for (const faPatch& fap : procMesh_.boundary())
412  {
413  for (const label patchEdgei : fap.edgeLabels())
414  {
415  const label globalEdgei = faEdgeProcAddr_[patchEdgei];
416 
417  const auto fnd = remapGlobal.cfind(globalEdgei);
418  if (fnd.found())
419  {
420  newEdgeProcAddr[edgei] = fnd.val();
421  ++edgei;
422  }
423  else
424  {
426  << "Failed to find edge " << globalEdgei
427  << " this indicates a programming error" << nl
428  << exit(FatalError);
429  }
430  }
431  }
432 
433  faEdgeProcAddr_ = std::move(newEdgeProcAddr);
434  }
435 }
436 
437 
438 void Foam::faMeshReconstructor::initPatch() const
439 {
440  serialPatchPtr_.reset
441  (
442  new primitivePatch
443  (
444  SubList<face>(singlePatchFaces_),
445  singlePatchPoints_
446  )
447  );
448 }
449 
450 
451 void Foam::faMeshReconstructor::createMesh()
452 {
453  const Time& runTime = procMesh_.mesh().time();
454 
455  // Time for non-parallel case (w/o functionObjects or libs)
456  serialRunTime_ = Time::New(runTime.globalPath().toAbsolute());
457 
458 
459  // Trivial polyMesh only containing points and faces.
460  // This is valid, provided we don't use it for anything complicated.
461 
462  serialVolMesh_.reset
463  (
464  new polyMesh
465  (
466  IOobject
467  (
468  procMesh_.mesh().name(), // Volume region name
469  procMesh_.mesh().facesInstance(),
470 
471  *serialRunTime_,
472 
475  false // not registered
476  ),
477  pointField(singlePatchPoints_), // copy
478  faceList(singlePatchFaces_), // copy
479  labelList(singlePatchFaces_.size(), Zero), // faceOwner
480  labelList(), // faceNeighbour
481  false // no syncPar!
482  )
483  );
484 
485  // A simple area-mesh with one-to-one mapping of faceLabels
486  serialAreaMesh_.reset
487  (
488  new faMesh
489  (
490  *serialVolMesh_,
491  identity(singlePatchFaces_.size()) // faceLabels
492  )
493  );
494 
495  auto& completeMesh = *serialAreaMesh_;
496 
497  // Add in non-processor boundary patches
498  PtrList<faPatch> completePatches(singlePatchEdgeLabels_.size());
499  forAll(completePatches, patchi)
500  {
501  const labelList& patchEdgeLabels = singlePatchEdgeLabels_[patchi];
502 
503  const faPatch& fap = procMesh_.boundary()[patchi];
504 
505  const label neiPolyPatchId = fap.ngbPolyPatchIndex();
506 
507  completePatches.set
508  (
509  patchi,
510  fap.clone
511  (
512  completeMesh.boundary(),
513  patchEdgeLabels,
514  patchi, // index
515  neiPolyPatchId
516  )
517  );
518  }
519 
520  // Serial mesh - no parallel communication
521 
522  const bool oldParRun = Pstream::parRun(false);
523 
524  completeMesh.addFaPatches(completePatches);
525 
526  Pstream::parRun(oldParRun); // Restore parallel state
527 }
528 
529 
530 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
531 
532 Foam::faMeshReconstructor::faMeshReconstructor
533 (
534  const faMesh& procMesh
535 )
536 :
537  procMesh_(procMesh)
538 {
539  if (!Pstream::parRun())
540  {
542  << "Can only be called in parallel!!" << nl
543  << exit(FatalError);
544  }
545 
546  // Require faceProcAddressing from finiteVolume decomposition
547  labelIOList fvFaceProcAddressing
548  (
549  IOobject
550  (
551  "faceProcAddressing",
552  procMesh_.mesh().facesInstance(),
553 
554  // Or search?
555  // procMesh_.time().findInstance
556  // (
557  // // Search for polyMesh face instance
558  // // mesh.facesInstance()
559  // procMesh_.mesh().meshDir(),
560  // "faceProcAddressing"
561  // ),
562 
564  procMesh_.mesh(), // The polyMesh db
567  false // not registered
568  )
569  );
570 
571  calcAddressing(fvFaceProcAddressing);
572 }
573 
574 
575 Foam::faMeshReconstructor::faMeshReconstructor
576 (
577  const faMesh& procMesh,
578  const labelUList& fvFaceProcAddressing
579 )
580 :
581  procMesh_(procMesh)
582 {
583  if (!Pstream::parRun())
584  {
586  << "Can only be called in parallel!!" << nl
587  << exit(FatalError);
588  }
589 
590  calcAddressing(fvFaceProcAddressing);
591 }
592 
593 
594 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
595 
597 {}
598 
599 
600 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
601 
603 {
604  serialAreaMesh_.reset(nullptr);
605  serialVolMesh_.reset(nullptr);
606  serialRunTime_.reset(nullptr);
607 }
608 
609 
611 {
612  if (!serialPatchPtr_)
613  {
614  initPatch();
615  }
616 
617  return *serialPatchPtr_;
618 }
619 
620 
622 {
623  if (!serialPatchPtr_)
624  {
625  initPatch();
626  }
627 
628  return *serialPatchPtr_;
629 }
630 
631 
633 {
634  if (!serialAreaMesh_)
635  {
636  const_cast<faMeshReconstructor&>(*this).createMesh();
637  }
638 
639  return *serialAreaMesh_;
640 }
641 
642 
644 {
645  writeAddressing(procMesh_.mesh().facesInstance());
646 }
647 
648 
650 {
651  // Write copies
652 
653  IOobject ioAddr
654  (
655  "procAddressing",
656  timeName,
658  procMesh_.mesh(), // The polyMesh
661  false // not registered
662  );
663 
664  // boundaryProcAddressing
665  ioAddr.rename("boundaryProcAddressing");
666  labelIOList(ioAddr, faBoundaryProcAddr_).write();
667 
668  // faceProcAddressing
669  ioAddr.rename("faceProcAddressing");
670  labelIOList(ioAddr, faFaceProcAddr_).write();
671 
672  // pointProcAddressing
673  ioAddr.rename("pointProcAddressing");
674  labelIOList(ioAddr, faPointProcAddr_).write();
675 
676  // edgeProcAddressing
677  ioAddr.rename("edgeProcAddressing");
678  labelIOList(ioAddr, faEdgeProcAddr_).write();
679 }
680 
681 
683 {
684  writeMesh(procMesh_.mesh().facesInstance());
685 }
686 
687 
689 {
690  const faMesh& fullMesh = this->mesh();
691 
692  const bool oldDistributed = fileHandler().distributed();
693  auto oldHandler = fileHandler(fileOperation::NewUncollated());
694  fileHandler().distributed(true);
695 
696  if (Pstream::master())
697  {
698  const bool oldParRun = Pstream::parRun(false);
699 
700  IOobject io(fullMesh.boundary());
701 
702  io.rename("faceLabels");
703  labelIOList(io, singlePatchFaceLabels_).write();
704 
705  fullMesh.boundary().write();
706 
707  Pstream::parRun(oldParRun);
708  }
709 
710  // Restore old settings
711  if (oldHandler)
712  {
713  fileHandler(std::move(oldHandler));
714  }
715  fileHandler().distributed(oldDistributed);
716 }
717 
718 
719 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::faMeshReconstructor::mesh
const faMesh & mesh() const
Serial equivalent faMesh.
Definition: faMeshReconstructor.C:632
Foam::faMeshReconstructor
A bare-bones reconstructor for finiteArea meshes when processor meshes are available (in parallel) bu...
Definition: faMeshReconstructor.H:64
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::faMesh::faceLabels
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:116
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
globalMeshData.H
Foam::faMeshReconstructor::~faMeshReconstructor
~faMeshReconstructor()
Destructor.
Definition: faMeshReconstructor.C:596
globalIndex.H
Foam::faMeshReconstructor::writeAddressing
void writeAddressing() const
Write proc addressing at the polyMesh faceInstances time.
Definition: faMeshReconstructor.C:643
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:321
Foam::fileOperation::NewUncollated
static autoPtr< fileOperation > NewUncollated()
Static construct the commonly used uncollatedFileOperation.
Definition: fileOperation.C:1474
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::combineReduce
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: PstreamCombineReduceOps.H:54
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1485
Foam::fileName::toAbsolute
fileName & toAbsolute()
Convert from relative to absolute.
Definition: fileName.C:377
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:44
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::globalMeshData::mergePoints
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
Definition: globalMeshData.C:2374
Foam::faMesh::nFaces
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:80
PstreamCombineReduceOps.H
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
nPatches
const label nPatches
Definition: printMeshSummary.H:30
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:515
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::faMeshReconstructor::patch
const primitivePatch & patch() const
Serial equivalent patch.
Definition: faMeshReconstructor.C:610
Foam::faMesh::nEdges
label nEdges() const noexcept
Number of local mesh edges.
Definition: faMeshI.H:62
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
Foam::faMesh::nInternalEdges
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:68
faMeshReconstructor.H
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::faMeshReconstructor::clearGeom
void clearGeom()
Definition: faMeshReconstructor.C:602
Foam::List::subList
SubList< T > subList
Declare type of subList.
Definition: List.H:112
Foam::faMesh::mesh
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::faMesh::boundary
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
Foam::Time::New
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Foam::IOobject::rename
virtual void rename(const word &newName)
Rename.
Definition: IOobject.H:506
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:508
edgeHashes.H
Foam::UPstream::commsTypes::nonBlocking
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::faBoundaryMesh::nNonProcessor
label nNonProcessor() const
The number of patches before the first processor patch.
Definition: faBoundaryMesh.C:188
Foam::UPstream::commsTypes::scheduled
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::faMesh::meshSubDir
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:471
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::IOList< label >
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::List::set
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:341
Foam::TimePaths::globalPath
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:82
Foam::primitivePatch
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Definition: primitivePatch.H:51
Foam::uindirectPrimitivePatch
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Definition: uindirectPrimitivePatch.H:49
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::faMesh::patch
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:122
Foam::faMeshReconstructor::debug
static int debug
Debug flag.
Definition: faMeshReconstructor.H:138
Foam::fileOperation::distributed
bool distributed() const noexcept
Distributed roots (parallel run)
Definition: fileOperation.H:239
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::faMeshReconstructor::writeMesh
void writeMesh() const
Write equivalent mesh information at the polyMesh faceInstances time.
Definition: faMeshReconstructor.C:682