reconstructParMesh.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2025 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  reconstructParMesh
29 
30 Group
31  grpParallelUtilities
32 
33 Description
34  Reconstructs a mesh using geometric information only.
35 
36  Writes point/face/cell procAddressing so afterwards reconstructPar can be
37  used to reconstruct fields.
38 
39 Usage
40  \b reconstructParMesh [OPTION]
41 
42  Options:
43  - \par -fullMatch
44  Does geometric matching on all boundary faces. Assumes no point
45  ordering
46 
47  - \par -procMatch
48  Assumes processor patches already in face order but not point order.
49  This is the pre v2106 default behaviour but might be removed if the new
50  topological method works well
51 
52  - \par -mergeTol <tol>
53  Specifies non-default merge tolerance (fraction of mesh bounding box)
54  for above options
55 
56  The default is to assume all processor boundaries are correctly ordered
57  (both faces and points) in which case no merge tolerance is needed.
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #include "argList.H"
62 #include "timeSelector.H"
63 
64 #include "IOobjectList.H"
65 #include "labelIOList.H"
66 #include "processorPolyPatch.H"
67 #include "mapAddedPolyMesh.H"
68 #include "polyMeshAdder.H"
69 #include "faceCoupleInfo.H"
70 #include "fvMeshAdder.H"
71 #include "polyTopoChange.H"
72 #include "topoSet.H"
73 #include "regionProperties.H"
74 #include "fvMeshTools.H"
75 
76 #include "faMeshReconstructor.H"
77 #include "ignoreFaPatch.H"
78 
79 using namespace Foam;
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
84 // usually meshes get written with limited precision (6 digits)
85 static const scalar defaultMergeTol = 1e-7;
86 
87 
88 // Determine which faces are coupled. Uses geometric merge distance.
89 // Looks either at all boundaryFaces (fullMatch) or only at the
90 // procBoundaries for proci. Assumes that masterMesh contains already merged
91 // all the processors < proci.
92 autoPtr<faceCoupleInfo> determineCoupledFaces
93 (
94  const bool fullMatch,
95  const label masterMeshProcStart,
96  const label masterMeshProcEnd,
97  const polyMesh& masterMesh,
98  const label meshToAddProcStart,
99  const label meshToAddProcEnd,
100  const polyMesh& meshToAdd,
101  const scalar mergeDist
102 )
103 {
104  if (fullMatch || masterMesh.nCells() == 0)
105  {
107  (
108  masterMesh,
109  meshToAdd,
110  mergeDist, // Absolute merging distance
111  true // Matching faces identical
112  );
113  }
114  else
115  {
116  // Pick up all patches on masterMesh ending in "toDDD" where DDD is
117  // the processor number proci.
118 
119  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
120 
121 
122  DynamicList<label> masterFaces
123  (
124  masterMesh.nFaces()
125  - masterMesh.nInternalFaces()
126  );
127 
128 
129  forAll(masterPatches, patchi)
130  {
131  const polyPatch& pp = masterPatches[patchi];
132 
133  if (isA<processorPolyPatch>(pp))
134  {
135  for
136  (
137  label proci=meshToAddProcStart;
138  proci<meshToAddProcEnd;
139  proci++
140  )
141  {
142  const string toProcString("to" + name(proci));
143  if (
144  pp.name().rfind(toProcString)
145  == (pp.name().size()-toProcString.size())
146  )
147  {
148  label meshFacei = pp.start();
149  forAll(pp, i)
150  {
151  masterFaces.append(meshFacei++);
152  }
153  break;
154  }
155  }
156 
157  }
158  }
159  masterFaces.shrink();
160 
161 
162  // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
163  // where DDD is the processor number proci and YYY is < proci.
164 
165  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
166 
167  DynamicList<label> addFaces
168  (
169  meshToAdd.nFaces()
170  - meshToAdd.nInternalFaces()
171  );
172 
173  forAll(addPatches, patchi)
174  {
175  const polyPatch& pp = addPatches[patchi];
176 
177  if (isA<processorPolyPatch>(pp))
178  {
179  bool isConnected = false;
180 
181  for
182  (
183  label mergedProci=masterMeshProcStart;
184  !isConnected && (mergedProci < masterMeshProcEnd);
185  mergedProci++
186  )
187  {
188  for
189  (
190  label proci = meshToAddProcStart;
191  proci < meshToAddProcEnd;
192  proci++
193  )
194  {
195  const word fromProcString
196  (
197  processorPolyPatch::newName(proci, mergedProci)
198  );
199 
200  if (pp.name() == fromProcString)
201  {
202  isConnected = true;
203  break;
204  }
205  }
206  }
207 
208  if (isConnected)
209  {
210  label meshFacei = pp.start();
211  forAll(pp, i)
212  {
213  addFaces.append(meshFacei++);
214  }
215  }
216  }
217  }
218  addFaces.shrink();
219 
221  (
222  masterMesh,
223  masterFaces,
224  meshToAdd,
225  addFaces,
226  mergeDist, // Absolute merging distance
227  true, // Matching faces identical?
228  false, // If perfect match are faces already ordered
229  // (e.g. processor patches)
230  false // are faces each on separate patch?
231  );
232  }
233 }
234 
235 
236 autoPtr<mapPolyMesh> mergeSharedPoints
237 (
238  const scalar mergeDist,
239  polyMesh& mesh,
240  labelListList& pointProcAddressing
241 )
242 {
243  // Find out which sets of points get merged and create a map from
244  // mesh point to unique point.
245  Map<label> pointToMaster
246  (
248  (
249  mesh,
250  mergeDist
251  )
252  );
253 
254  Info<< "mergeSharedPoints : detected " << pointToMaster.size()
255  << " points that are to be merged." << endl;
256 
257  if (returnReduceAnd(pointToMaster.empty()))
258  {
259  return nullptr;
260  }
261 
262  polyTopoChange meshMod(mesh);
263 
264  fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
265 
266  // Change the mesh (no inflation). Note: parallel comms allowed.
267  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
268 
269  // Update fields. No inflation, parallel sync.
270  mesh.updateMesh(map());
271 
272  // pointProcAddressing give indices into the master mesh so adapt them
273  // for changed point numbering.
274 
275  // Adapt constructMaps for merged points.
276  forAll(pointProcAddressing, proci)
277  {
278  labelList& constructMap = pointProcAddressing[proci];
279 
280  forAll(constructMap, i)
281  {
282  label oldPointi = constructMap[i];
283 
284  // New label of point after changeMesh.
285  label newPointi = map().reversePointMap()[oldPointi];
286 
287  if (newPointi < -1)
288  {
289  constructMap[i] = -newPointi-2;
290  }
291  else if (newPointi >= 0)
292  {
293  constructMap[i] = newPointi;
294  }
295  else
296  {
298  << "Problem. oldPointi:" << oldPointi
299  << " newPointi:" << newPointi << abort(FatalError);
300  }
301  }
302  }
303 
304  return map;
305 }
306 
307 
308 boundBox procBounds
309 (
310  const PtrList<Time>& databases,
311  const word& regionName
312 )
313 {
314  boundBox bb;
315 
316  for (const Time& procDb : databases)
317  {
318  fileName pointsInstance
319  (
320  procDb.findInstance(polyMesh::meshDir(regionName), "points")
321  );
322 
323  if (pointsInstance != procDb.timeName())
324  {
326  << "Your time was specified as " << procDb.timeName()
327  << " but there is no polyMesh/points in that time." << nl
328  << "(points file at " << pointsInstance << ')' << nl
329  << "Please rerun with the correct time specified"
330  << " (through the -constant, -time or -latestTime "
331  << "(at your option)."
332  << endl << exit(FatalError);
333  }
334 
335  Info<< "Reading points from "
336  << procDb.caseName()
337  << " for time = " << procDb.timeName()
338  << nl << endl;
339 
341  (
342  IOobject
343  (
344  "points",
345  pointsInstance,
347  procDb,
351  )
352  );
353 
354  bb.add(points);
355  }
356 
357  return bb;
358 }
359 
360 
361 void writeDistribution
362 (
363  Time& runTime,
364  const fvMesh& masterMesh,
365  const labelListList& cellProcAddressing
366 )
367 {
368  // Write the decomposition as labelList for use with 'manual'
369  // decomposition method.
370  labelIOList cellDecomposition
371  (
372  IOobject
373  (
374  "cellDecomposition",
375  masterMesh.facesInstance(),
376  masterMesh,
380  ),
381  masterMesh.nCells()
382  );
383 
384  forAll(cellProcAddressing, proci)
385  {
386  const labelList& pCells = cellProcAddressing[proci];
387  labelUIndList(cellDecomposition, pCells) = proci;
388  }
389 
390  cellDecomposition.write();
391 
392  Info<< nl << "Wrote decomposition to "
393  << cellDecomposition.objectRelPath()
394  << " for use in manual decomposition." << endl;
395 
396  // Write as volScalarField for postprocessing.
397  // Change time to 0 if was 'constant'
398  {
399  const scalar oldTime = runTime.value();
400  const label oldIndex = runTime.timeIndex();
401  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
402  {
403  runTime.setTime(0, oldIndex+1);
404  }
405 
406  volScalarField cellDist
407  (
408  IOobject
409  (
410  "cellDist",
411  runTime.timeName(),
412  masterMesh,
416  ),
417  masterMesh,
420  );
421 
422  forAll(cellDecomposition, celli)
423  {
424  cellDist[celli] = cellDecomposition[celli];
425  }
426 
427  cellDist.correctBoundaryConditions();
428  cellDist.write();
429 
430  Info<< nl << "Wrote decomposition to "
431  << cellDist.objectRelPath()
432  << " (volScalarField) for visualization."
433  << endl;
434 
435  // Restore time
436  runTime.setTime(oldTime, oldIndex);
437  }
438 }
439 
440 
441 void writeMesh
442 (
443  const fvMesh& mesh,
444  const labelListList& cellProcAddressing
445 )
446 {
447  const fileName outputDir
448  (
449  mesh.time().path()
450  / mesh.time().timeName()
451  / mesh.regionName()
453  );
454 
455  Info<< nl << "Writing merged mesh to "
456  << mesh.time().relativePath(outputDir) << nl << endl;
457 
458  if (!mesh.write())
459  {
461  << "Failed writing polyMesh."
462  << exit(FatalError);
463  }
465 }
466 
467 
468 void writeMaps
469 (
470  const label masterInternalFaces,
471  const labelUList& masterOwner,
472  const polyMesh& procMesh,
473  const labelUList& cellProcAddressing,
474  const labelUList& faceProcAddressing,
475  const labelUList& pointProcAddressing,
476  const labelUList& boundProcAddressing
477 )
478 {
479  const fileName outputDir
480  (
481  procMesh.time().caseName()
482  / procMesh.facesInstance()
483  / procMesh.regionName()
485  );
486 
487  IOobject ioAddr
488  (
489  "procAddressing",
490  procMesh.facesInstance(),
492  procMesh,
496  );
497 
498 
499  Info<< "Writing addressing : " << outputDir << nl;
500 
501  // From processor point to reconstructed mesh point
502 
503  Info<< " pointProcAddressing" << endl;
504  ioAddr.rename("pointProcAddressing");
505  IOList<label>::writeContents(ioAddr, pointProcAddressing);
506 
507  // From processor face to reconstructed mesh face
508  Info<< " faceProcAddressing" << endl;
509  ioAddr.rename("faceProcAddressing");
510  labelIOList faceProcAddr(ioAddr, faceProcAddressing);
511 
512  // Now add turning index to faceProcAddressing.
513  // See reconstructPar for meaning of turning index.
514  forAll(faceProcAddr, procFacei)
515  {
516  const label masterFacei = faceProcAddr[procFacei];
517 
518  if
519  (
520  !procMesh.isInternalFace(procFacei)
521  && masterFacei < masterInternalFaces
522  )
523  {
524  // proc face is now external but used to be internal face.
525  // Check if we have owner or neighbour.
526 
527  label procOwn = procMesh.faceOwner()[procFacei];
528  label masterOwn = masterOwner[masterFacei];
529 
530  if (cellProcAddressing[procOwn] == masterOwn)
531  {
532  // No turning. Offset by 1.
533  faceProcAddr[procFacei]++;
534  }
535  else
536  {
537  // Turned face.
538  faceProcAddr[procFacei] = -1 - faceProcAddr[procFacei];
539  }
540  }
541  else
542  {
543  // No turning. Offset by 1.
544  faceProcAddr[procFacei]++;
545  }
546  }
547 
548  faceProcAddr.write();
549 
550 
551  // From processor cell to reconstructed mesh cell
552  Info<< " cellProcAddressing" << endl;
553  ioAddr.rename("cellProcAddressing");
554  IOList<label>::writeContents(ioAddr, cellProcAddressing);
555 
556 
557  // From processor patch to reconstructed mesh patch
558  Info<< " boundaryProcAddressing" << endl;
559  ioAddr.rename("boundaryProcAddressing");
560  IOList<label>::writeContents(ioAddr, boundProcAddressing);
561 
562  Info<< endl;
563 }
564 
565 
566 void printWarning()
567 {
568  Info<<
569 "Merge individual processor meshes back into one master mesh.\n"
570 "Use if the original master mesh has been deleted or the processor meshes\n"
571 "have been modified (topology change).\n"
572 "This tool will write the resulting mesh to a new time step and construct\n"
573 "xxxxProcAddressing files in the processor meshes so reconstructPar can be\n"
574 "used to regenerate the fields on the master mesh.\n\n"
575 "Not well tested & use at your own risk!\n\n";
576 }
577 
578 
579 // Used to determine the correspondence between the edges.
580 // See faMeshReconstructor::calcAddressing
581 void determineFaEdgeMapping
582 (
583  const uindirectPrimitivePatch& onePatch,// reconstructed faMesh patch
584  const PtrList<faMesh>& procFaMeshes, // individual faMeshes
585  const labelListList& pointProcAddressing, // procPolyMesh to reconstructed
586 
587  labelListList& faEdgeProcAddressing
588 )
589 {
590  // Determines ordering from faMesh to onePatch (= recsontructed faMesh)
591 
592  // From two polyMesh points to patch-edge label
593  EdgeMap<label> pointsToOnePatchEdge(onePatch.nEdges());
594  {
595  const edgeList& edges = onePatch.edges();
596  const labelList& mp = onePatch.meshPoints();
597 
598  forAll(edges, edgei)
599  {
600  const edge meshE(mp, edges[edgei]);
601  pointsToOnePatchEdge.insert(meshE, edgei);
602  }
603  }
604 
605  // Now we can create the mapping from patch-edge on processor
606  // to patch-edge on onePatch
607 
608  faEdgeProcAddressing.resize_nocopy(procFaMeshes.size());
609  forAll(procFaMeshes, proci)
610  {
611  const auto& procPatch = procFaMeshes[proci].patch();
612  const auto& edges = procPatch.edges();
613  const auto& mp = procPatch.meshPoints();
614  const auto& ppAddressing = pointProcAddressing[proci];
615 
616  labelList& edgeProcAddr = faEdgeProcAddressing[proci];
617  edgeProcAddr.resize_nocopy(edges.size());
618 
619  label edgei = 0;
620  for
621  (
622  ;
623  edgei < procPatch.nEdges(); //procPatch.nInternalEdges();
624  edgei++
625  )
626  {
627  const edge meshE(mp, edges[edgei]);
628  const edge onePatchE(ppAddressing, meshE);
629 
630  edgeProcAddr[edgei] = pointsToOnePatchEdge[onePatchE];
631  }
632  }
633 }
634 
635 
636 void sortFaEdgeMapping
637 (
638  const uindirectPrimitivePatch& onePatch,// reconstructed faMesh patch
639  const PtrList<faMesh>& procFaMeshes, // individual faMeshes
640  const labelListList& pointProcAddressing, // procPolyMesh to reconstructed
641 
642  labelListList& faEdgeProcAddressing, // per proc the map to the master
643  labelListList& singlePatchEdgeLabels // per patch the map to the master
644 )
645 {
646  // From faMeshReconstructor.C - edge shuffling on patches
647 
648  Map<label> remapGlobal(2*onePatch.nEdges());
649  for (label edgei = 0; edgei < onePatch.nInternalEdges(); ++edgei)
650  {
651  remapGlobal.insert(edgei, remapGlobal.size());
652  }
653 
654 
655  const faMesh& procMesh0 = procFaMeshes[0];
656  const label nNonProcessor = procMesh0.boundary().nNonProcessor();
657 
658  singlePatchEdgeLabels.resize_nocopy(nNonProcessor);
659 
660  forAll(singlePatchEdgeLabels, patchi)
661  {
662  if (isA<ignoreFaPatch>(procMesh0.boundary()[patchi]))
663  {
664  // These are not real edges
665  continue;
666  }
667 
668  forAll(procFaMeshes, proci)
669  {
670  const faMesh& procMesh = procFaMeshes[proci];
671  const faPatch& fap = procMesh.boundary()[patchi];
672 
673  labelList patchEdgeLabels(fap.edgeLabels());
674 
675  // Renumber from local edges to global edges (natural order)
676  for (label& edgeId : patchEdgeLabels)
677  {
678  edgeId = faEdgeProcAddressing[proci][edgeId];
679  }
680 
681  // Combine from all processors
682  singlePatchEdgeLabels[patchi].append(patchEdgeLabels);
683  }
684 
685  // Sorted order will be the original non-decomposed order
686  Foam::sort(singlePatchEdgeLabels[patchi]);
687 
688  for (const label sortedEdgei : singlePatchEdgeLabels[patchi])
689  {
690  remapGlobal.insert(sortedEdgei, remapGlobal.size());
691  }
692  }
693 
694  {
695  // Use the map to rewrite the local faEdgeProcAddressing
696 
697  labelListList newEdgeProcAddr(faEdgeProcAddressing);
698 
699  forAll(procFaMeshes, proci)
700  {
701  const faMesh& procMesh = procFaMeshes[proci];
702 
703  label edgei = procMesh.nInternalEdges();
704 
705  for (const faPatch& fap : procMesh.boundary())
706  {
707  for (const label patchEdgei : fap.edgeLabels())
708  {
709  const label globalEdgei =
710  faEdgeProcAddressing[proci][patchEdgei];
711 
712  const auto fnd = remapGlobal.cfind(globalEdgei);
713  if (fnd.good())
714  {
715  newEdgeProcAddr[proci][edgei] = fnd.val();
716  ++edgei;
717  }
718  else
719  {
721  << "Failed to find edge " << globalEdgei
722  << " this indicates a programming error" << nl
723  << exit(FatalError);
724  }
725  }
726  }
727  }
728  faEdgeProcAddressing = std::move(newEdgeProcAddr);
729  }
730 }
731 
732 
733 int main(int argc, char *argv[])
734 {
736  (
737  "Reconstruct a mesh using geometric/topological information only"
738  );
739 
740  // Enable -constant ... if someone really wants it
741  // Enable -withZero to prevent accidentally trashing the initial fields
742  timeSelector::addOptions(true, true); // constant(true), zero(true)
743 
745 
748  (
749  "addressing-only",
750  "Create procAddressing only without overwriting the mesh"
751  );
753  (
754  "mergeTol",
755  "scalar",
756  "The merge distance relative to the bounding box size (default 1e-7)"
757  );
759  (
760  "fullMatch",
761  "Do (slower) geometric matching on all boundary faces"
762  );
764  (
765  "procMatch",
766  "Do matching on processor faces only"
767  );
769  (
770  "cellDist",
771  "Write cell distribution as a labelList - for use with 'manual' "
772  "decomposition method or as a volScalarField for post-processing."
773  );
775  (
776  "no-finite-area",
777  "Suppress finiteArea mesh reconstruction",
778  true // Advanced option
779  );
780 
781  #include "addAllRegionOptions.H"
782 
783  #include "setRootCase.H"
784  #include "createTime.H"
785 
786  printWarning();
787 
788  const bool fullMatch = args.found("fullMatch");
789  const bool procMatch = args.found("procMatch");
790  const bool writeCellDist = args.found("cellDist");
791  bool doFiniteArea = !args.found("no-finite-area");
792  const bool writeAddrOnly = args.found("addressing-only");
793 
794  const scalar mergeTol =
795  args.getOrDefault<scalar>("mergeTol", defaultMergeTol);
796 
797  if (fullMatch)
798  {
799  Info<< "Use geometric matching on all boundary faces." << nl << endl;
800  }
801  else if (procMatch)
802  {
803  Info<< "Use geometric matching on correct procBoundaries only." << nl
804  << "This assumes a correct decomposition." << endl;
805  }
806  else
807  {
808  Info<< "Merge assuming correct, fully matched procBoundaries." << nl
809  << endl;
810  }
811 
812  if (fullMatch || procMatch)
813  {
814  const scalar writeTol =
815  std::pow(scalar(10), -scalar(IOstream::defaultPrecision()));
816 
817  Info<< "Merge tolerance : " << mergeTol << nl
818  << "Write tolerance : " << writeTol << endl;
819 
820  if
821  (
823  && mergeTol < writeTol
824  )
825  {
827  << "Your current settings specify ASCII writing with "
828  << IOstream::defaultPrecision() << " digits precision." << endl
829  << "Your merging tolerance (" << mergeTol << ")"
830  << " is finer than this." << endl
831  << "Please change your writeFormat to binary"
832  << " or increase the writePrecision" << endl
833  << "or adjust the merge tolerance (-mergeTol)."
834  << exit(FatalError);
835  }
836  }
837 
838  // Get region names
839  #include "getAllRegionOptions.H"
840 
841  // Determine the processor count
842  label nProcs{0};
843 
844  if (regionNames.empty())
845  {
847  << "No regions specified or detected."
848  << exit(FatalError);
849  }
850  else if (regionNames[0] == polyMesh::defaultRegion)
851  {
852  nProcs = fileHandler().nProcs(args.path());
853  }
854  else
855  {
856  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
857 
858  if (regionNames.size() == 1)
859  {
860  Info<< "Using region: " << regionNames[0] << nl << endl;
861  }
862  }
863 
864  if (!nProcs)
865  {
867  << "No processor* directories found"
868  << exit(FatalError);
869  }
870 
871  // Read all time databases
872  PtrList<Time> databases(nProcs);
873 
874  Info<< "Found " << nProcs << " processor directories" << endl;
875  forAll(databases, proci)
876  {
877  Info<< " Reading database "
878  << args.caseName()/("processor" + Foam::name(proci))
879  << endl;
880 
881  databases.set
882  (
883  proci,
884  new Time
885  (
887  args.rootPath(),
888  args.caseName()/("processor" + Foam::name(proci)),
890  args.allowLibs()
891  )
892  );
893  }
894  Info<< endl;
895 
896  // Use the times list from the master processor
897  // and select a subset based on the command-line options
899  (
900  databases[0].times(),
901  args
902  );
903 
904  // Loop over all times
905  forAll(timeDirs, timei)
906  {
907  // Set time for global database
908  runTime.setTime(timeDirs[timei], timei);
909 
910  // Set time for all databases
911  forAll(databases, proci)
912  {
913  databases[proci].setTime(timeDirs[timei], timei);
914  }
915 
916  Info<< "Time = " << runTime.timeName() << endl;
917 
918  // Check for any mesh changes
919  label nMeshChanged = 0;
920  boolList hasRegionMesh(regionNames.size(), false);
921  forAll(regionNames, regioni)
922  {
923  const word& regionName = regionNames[regioni];
924 
925  IOobject facesIO
926  (
927  "faces",
928  databases[0].timeName(),
930  databases[0],
933  );
934 
935  // Problem: faceCompactIOList recognises both 'faceList' and
936  // 'faceCompactList' so we should be lenient when doing
937  // typeHeaderOk
938 
939  const bool ok = facesIO.typeHeaderOk<faceCompactIOList>(false);
940 
941  if (ok)
942  {
943  hasRegionMesh[regioni] = true;
944  ++nMeshChanged;
945  }
946  }
947 
948  // Check for any mesh changes
949  if (!nMeshChanged)
950  {
951  Info<< "No meshes" << nl << endl;
952  continue;
953  }
954 
955  Info<< endl;
956 
957  forAll(regionNames, regioni)
958  {
959  const word& regionName = regionNames[regioni];
961 
962  if (!hasRegionMesh[regioni])
963  {
964  Info<< "region=" << regionName << " (no mesh)" << nl << endl;
965  continue;
966  }
967 
968  if (regionNames.size() > 1)
969  {
970  Info<< "region=" << regionName << nl;
971  }
972 
973 
974  // Addressing from processor to reconstructed case
975  labelListList cellProcAddressing(nProcs);
976  labelListList faceProcAddressing(nProcs);
977  labelListList pointProcAddressing(nProcs);
978  labelListList boundProcAddressing(nProcs);
979 
980 
981  // Internal faces on the final reconstructed mesh
982  label masterInternalFaces;
983 
984  // Owner addressing on the final reconstructed mesh
985  labelList masterOwner;
986 
987  if (procMatch)
988  {
989  // Read point on individual processors to determine
990  // merge tolerance
991  // (otherwise single cell domains might give problems)
992 
993  const boundBox bb = procBounds(databases, regionDir);
994  const scalar mergeDist = mergeTol*bb.mag();
995 
996  Info<< "Overall mesh bounding box : " << bb << nl
997  << "Relative tolerance : " << mergeTol << nl
998  << "Absolute matching distance : " << mergeDist << nl
999  << endl;
1000 
1001 
1002  // Construct empty mesh.
1003  PtrList<fvMesh> masterMesh(nProcs);
1004 
1005  for (label proci=0; proci<nProcs; proci++)
1006  {
1007  masterMesh.set
1008  (
1009  proci,
1010  new fvMesh
1011  (
1012  IOobject
1013  (
1014  regionName,
1015  runTime.timeName(),
1016  runTime,
1018  ),
1019  Zero
1020  )
1021  );
1022 
1023  fvMesh meshToAdd
1024  (
1025  IOobject
1026  (
1027  regionName,
1028  databases[proci].timeName(),
1029  databases[proci]
1030  )
1031  );
1032 
1033  // Initialize its addressing
1034  cellProcAddressing[proci] = identity(meshToAdd.nCells());
1035  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
1036  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
1037  boundProcAddressing[proci] =
1038  identity(meshToAdd.boundaryMesh().size());
1039 
1040  // Find geometrically shared points/faces.
1041  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
1042  (
1043  fullMatch,
1044  proci,
1045  proci,
1046  masterMesh[proci],
1047  proci,
1048  proci,
1049  meshToAdd,
1050  mergeDist
1051  );
1052 
1053  // Add elements to mesh
1055  (
1056  masterMesh[proci],
1057  meshToAdd,
1058  couples()
1059  );
1060 
1061  // Added processor
1062  renumber(map().addedCellMap(), cellProcAddressing[proci]);
1063  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
1064  renumber(map().addedPointMap(), pointProcAddressing[proci]);
1065  renumber(map().addedPatchMap(), boundProcAddressing[proci]);
1066  }
1067  for (label step=2; step<nProcs*2; step*=2)
1068  {
1069  for (label proci=0; proci<nProcs; proci+=step)
1070  {
1071  label next = proci + step/2;
1072  if(next >= nProcs)
1073  {
1074  continue;
1075  }
1076 
1077  Info<< "Merging mesh " << proci << " with "
1078  << next << endl;
1079 
1080  // Find geometrically shared points/faces.
1081  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
1082  (
1083  fullMatch,
1084  proci,
1085  next,
1086  masterMesh[proci],
1087  next,
1088  proci+step,
1089  masterMesh[next],
1090  mergeDist
1091  );
1092 
1093  // Add elements to mesh
1095  (
1096  masterMesh[proci],
1097  masterMesh[next],
1098  couples()
1099  );
1100 
1101  // Processors that were already in masterMesh
1102  for (label mergedI=proci; mergedI<next; mergedI++)
1103  {
1104  renumber
1105  (
1106  map().oldCellMap(),
1107  cellProcAddressing[mergedI]
1108  );
1109 
1110  renumber
1111  (
1112  map().oldFaceMap(),
1113  faceProcAddressing[mergedI]
1114  );
1115 
1116  renumber
1117  (
1118  map().oldPointMap(),
1119  pointProcAddressing[mergedI]
1120  );
1121 
1122  // Note: boundary is special since can contain -1.
1123  renumber
1124  (
1125  map().oldPatchMap(),
1126  boundProcAddressing[mergedI]
1127  );
1128  }
1129 
1130  // Added processor
1131  for
1132  (
1133  label addedI=next;
1134  addedI < Foam::min(proci+step, nProcs);
1135  addedI++
1136  )
1137  {
1138  renumber
1139  (
1140  map().addedCellMap(),
1141  cellProcAddressing[addedI]
1142  );
1143 
1144  renumber
1145  (
1146  map().addedFaceMap(),
1147  faceProcAddressing[addedI]
1148  );
1149 
1150  renumber
1151  (
1152  map().addedPointMap(),
1153  pointProcAddressing[addedI]
1154  );
1155 
1156  renumber
1157  (
1158  map().addedPatchMap(),
1159  boundProcAddressing[addedI]
1160  );
1161  }
1162 
1163  masterMesh.set(next, nullptr);
1164  }
1165  }
1166 
1167  for (label proci=0; proci<nProcs; proci++)
1168  {
1169  Info<< "Reading mesh to add from "
1170  << databases[proci].caseName()
1171  << " for time = " << databases[proci].timeName()
1172  << nl << nl << endl;
1173  }
1174 
1175  // See if any points on the mastermesh have become connected
1176  // because of connections through processor meshes.
1177  mergeSharedPoints(mergeDist,masterMesh[0],pointProcAddressing);
1178 
1179  // Save some properties on the reconstructed mesh
1180  masterInternalFaces = masterMesh[0].nInternalFaces();
1181  masterOwner = masterMesh[0].faceOwner();
1182 
1183 
1184  if (writeAddrOnly)
1185  {
1186  Info<< nl
1187  << "Disabled writing of merged mesh (-addressing-only)"
1188  << nl << nl;
1189  }
1190  else
1191  {
1192  // Write reconstructed mesh
1193  writeMesh(masterMesh[0], cellProcAddressing);
1194  if (writeCellDist)
1195  {
1196  writeDistribution
1197  (
1198  runTime,
1199  masterMesh[0],
1200  cellProcAddressing
1201  );
1202  }
1203  }
1204  }
1205  else
1206  {
1207  // Load all meshes
1208  PtrList<fvMesh> fvMeshes(nProcs);
1209  for (label proci=0; proci<nProcs; proci++)
1210  {
1211  fvMeshes.set
1212  (
1213  proci,
1214  new fvMesh
1215  (
1216  IOobject
1217  (
1218  regionName,
1219  databases[proci].timeName(),
1220  databases[proci]
1221  )
1222  )
1223  );
1224  }
1225 
1226  // Construct pointers to meshes
1227  UPtrList<polyMesh> meshes(fvMeshes.size());
1228  forAll(fvMeshes, proci)
1229  {
1230  meshes.set(proci, &fvMeshes[proci]);
1231  }
1232 
1233  // Get pairs of patches to stitch. These pairs have to
1234  // - have ordered, opposite faces (so one to one correspondence)
1235  List<DynamicList<label>> localPatch;
1236  List<DynamicList<label>> remoteProc;
1237  List<DynamicList<label>> remotePatch;
1238  const label nGlobalPatches = polyMeshAdder::procPatchPairs
1239  (
1240  meshes,
1241  localPatch,
1242  remoteProc,
1243  remotePatch
1244  );
1245 
1246  // Collect matching boundary faces on patches-to-stitch
1247  labelListList localBoundaryFace;
1248  labelListList remoteFaceProc;
1249  labelListList remoteBoundaryFace;
1251  (
1252  meshes,
1253  localPatch,
1254  remoteProc,
1255  remotePatch,
1256  localBoundaryFace,
1257  remoteFaceProc,
1258  remoteBoundaryFace
1259  );
1260 
1261  // All matched faces assumed to have vertex0 matched
1262  labelListList remoteFaceStart(meshes.size());
1263  forAll(meshes, proci)
1264  {
1265  const labelList& procFaces = localBoundaryFace[proci];
1266  remoteFaceStart[proci].setSize(procFaces.size(), 0);
1267  }
1268 
1269 
1270 
1271  labelListList patchMap(meshes.size());
1272  labelListList pointZoneMap(meshes.size());
1273  labelListList faceZoneMap(meshes.size());
1274  labelListList cellZoneMap(meshes.size());
1275  forAll(meshes, proci)
1276  {
1277  const polyMesh& mesh = meshes[proci];
1278  patchMap[proci] = identity(mesh.boundaryMesh().size());
1279 
1280  // Remove excess patches
1281  patchMap[proci].setSize(nGlobalPatches);
1282 
1283  pointZoneMap[proci] = identity(mesh.pointZones().size());
1284  faceZoneMap[proci] = identity(mesh.faceZones().size());
1285  cellZoneMap[proci] = identity(mesh.cellZones().size());
1286  }
1287 
1288 
1289  refPtr<fvMesh> masterMeshPtr;
1290  {
1291  // Do in-place addition on proc0.
1292 
1293  const labelList oldFaceOwner(fvMeshes[0].faceOwner());
1294 
1296  (
1297  0, // index of mesh to modify (== mesh_)
1298  fvMeshes,
1299  oldFaceOwner,
1300 
1301  // Coupling info
1302  localBoundaryFace,
1303  remoteFaceProc,
1304  remoteBoundaryFace,
1305 
1306  boundProcAddressing,
1307  cellProcAddressing,
1308  faceProcAddressing,
1309  pointProcAddressing
1310  );
1311 
1312  // Remove zero-faces processor patches
1313  const polyBoundaryMesh& pbm = fvMeshes[0].boundaryMesh();
1314  labelList oldToNew(pbm.size(), -1);
1315  label newi = 0;
1316  // Non processor patches first
1317  forAll(pbm, patchi)
1318  {
1319  const auto& pp = pbm[patchi];
1320  if (!isA<processorPolyPatch>(pp) || pp.size())
1321  {
1322  oldToNew[patchi] = newi++;
1323  }
1324  }
1325  const label nNonProcPatches = newi;
1326 
1327  // Move all deletable patches to the end
1328  forAll(oldToNew, patchi)
1329  {
1330  if (oldToNew[patchi] == -1)
1331  {
1332  oldToNew[patchi] = newi++;
1333  }
1334  }
1336  (
1337  fvMeshes[0],
1338  oldToNew,
1339  nNonProcPatches,
1340  false
1341  );
1342 
1343  masterMeshPtr.cref(fvMeshes[0]);
1344  }
1345 
1346 
1347  const fvMesh& masterMesh = masterMeshPtr();
1348 
1349  // Number of internal faces on the final reconstructed mesh
1350  masterInternalFaces = masterMesh.nInternalFaces();
1351 
1352  // Owner addressing on the final reconstructed mesh
1353  masterOwner = masterMesh.faceOwner();
1354 
1355 
1356  // Write reconstructed mesh
1357  // Override:
1358  // - caseName
1359  // - processorCase flag
1360  // so the resulting mesh goes to the correct location (even with
1361  // collated). The better way of solving this is to construct
1362  // (zero) mesh on the undecomposed runTime.
1363 
1364  if (writeAddrOnly)
1365  {
1366  Info<< nl
1367  << "Disabled writing of merged mesh (-addressing-only)"
1368  << nl << nl;
1369  }
1370  else
1371  {
1372  Time& masterTime = const_cast<Time&>(masterMesh.time());
1373 
1374  const word oldCaseName = masterTime.caseName();
1375  masterTime.caseName() = runTime.caseName();
1376  const bool oldProcCase(masterTime.processorCase(false));
1377 
1378  // Write reconstructed mesh
1379  writeMesh(masterMesh, cellProcAddressing);
1380  if (writeCellDist)
1381  {
1382  writeDistribution
1383  (
1384  runTime,
1385  masterMesh,
1386  cellProcAddressing
1387  );
1388  }
1389  masterTime.caseName() = oldCaseName;
1390  masterTime.processorCase(oldProcCase);
1391  }
1392  }
1393 
1394 
1395  // Write the addressing
1396 
1397  Info<< "Reconstructing addressing from processor meshes"
1398  << " to the newly reconstructed mesh" << nl << endl;
1399 
1400 
1401  // Read finite-area
1402  PtrList<faMesh> procFaMeshes(databases.size());
1403  PtrList<polyMesh> procMeshes(databases.size());
1404 
1405  forAll(databases, proci)
1406  {
1407  Info<< "Processor " << proci << nl
1408  << "Read processor mesh: "
1409  << (databases[proci].caseName()/regionDir) << endl;
1410 
1411  procMeshes.set
1412  (
1413  proci,
1414  new polyMesh
1415  (
1416  IOobject
1417  (
1418  regionName,
1419  databases[proci].timeName(),
1420  databases[proci]
1421  )
1422  )
1423  );
1424  const polyMesh& procMesh = procMeshes[proci];
1425 
1426  writeMaps
1427  (
1428  masterInternalFaces,
1429  masterOwner,
1430  procMesh,
1431  cellProcAddressing[proci],
1432  faceProcAddressing[proci],
1433  pointProcAddressing[proci],
1434  boundProcAddressing[proci]
1435  );
1436 
1437 
1438  if (doFiniteArea)
1439  {
1440  const word boundaryInst =
1441  procMesh.time().findInstance
1442  (
1443  procMesh.meshDir(),
1444  "boundary"
1445  );
1446 
1447  IOobject io
1448  (
1449  "faBoundary",
1450  boundaryInst,
1451  faMesh::meshDir(procMesh, word::null),
1452  procMesh.time(),
1456  );
1457 
1458  if (io.typeHeaderOk<faBoundaryMesh>(true))
1459  {
1460  // Always based on the volume decomposition!
1461  procFaMeshes.set(proci, new faMesh(procMesh));
1462  }
1463  }
1464  }
1465 
1466 
1467  // Finite-area mapping
1468  doFiniteArea = false;
1469  forAll(procFaMeshes, proci)
1470  {
1471  if (procFaMeshes.set(proci))
1472  {
1473  doFiniteArea = true;
1474  }
1475  }
1476 
1477  if (doFiniteArea)
1478  {
1479  // Addressing from processor to reconstructed case
1480  labelListList faFaceProcAddressing(nProcs);
1481  labelListList faEdgeProcAddressing(nProcs);
1482  labelListList faPointProcAddressing(nProcs);
1483  labelListList faBoundProcAddressing(nProcs);
1484 
1485 
1486  // boundProcAddressing
1487  // ~~~~~~~~~~~~~~~~~~~
1488 
1489  forAll(procFaMeshes, proci)
1490  {
1491  const auto& procMesh = procFaMeshes[proci];
1492  const auto& bm = procMesh.boundary();
1493 
1494  faBoundProcAddressing[proci] = identity(bm.size());
1495  // Mark processor patches
1496  for
1497  (
1498  label patchi = bm.nNonProcessor();
1499  patchi < bm.size();
1500  ++patchi
1501  )
1502  {
1503  faBoundProcAddressing[proci][patchi] = -1;
1504  }
1505  }
1506 
1507 
1508  // Re-read reconstructed polyMesh. Note: could probably be avoided
1509  // by merging into loops above.
1510  const polyMesh masterMesh
1511  (
1512  IOobject
1513  (
1514  regionName,
1515  runTime.timeName(),
1516  runTime
1517  ),
1518  true
1519  );
1520 
1521 
1522  // faceProcAddressing
1523  // ~~~~~~~~~~~~~~~~~~
1524 
1525  DynamicList<label> masterFaceLabels;
1526  label nPatchFaces = 0;
1527  forAll(procFaMeshes, proci)
1528  {
1529  const auto& procMesh = procFaMeshes[proci];
1530  const auto& procPolyFaces = procMesh.faceLabels();
1531  const auto& fpa = faceProcAddressing[proci];
1532 
1533  labelList& faceAddr = faFaceProcAddressing[proci];
1534  faceAddr.resize_nocopy(procPolyFaces.size());
1535 
1536  // Map to masterPolyMesh faces
1537  forAll(procPolyFaces, i)
1538  {
1539  const label facei = procPolyFaces[i];
1540  masterFaceLabels.append(fpa[facei]);
1541  faceAddr[i] = nPatchFaces++;
1542  }
1543  }
1544 
1545 
1546  // faMesh itself
1547  // ~~~~~~~~~~~~~
1548 
1549  // Set up to read-if-present
1550  IOobject io(masterMesh);
1552 
1553  // Construct without patches
1554  faMesh masterFaMesh
1555  (
1556  masterMesh,
1557  std::move(masterFaceLabels),
1558  io
1559  );
1560 
1561  const uindirectPrimitivePatch& masterPatch =
1562  masterFaMesh.patch();
1563 
1564 
1565  // pointProcAddressing
1566  // ~~~~~~~~~~~~~~~~~~~
1567 
1568  const auto& mpm = masterPatch.meshPointMap();
1569  forAll(procFaMeshes, proci)
1570  {
1571  const auto& procPatch = procFaMeshes[proci].patch();
1572  const auto& mp = procPatch.meshPoints();
1573 
1574  labelList& pointAddr = faPointProcAddressing[proci];
1575  pointAddr.resize_nocopy(mp.size());
1576 
1577  forAll(mp, i)
1578  {
1579  pointAddr[i] = mpm[pointProcAddressing[proci][mp[i]]];
1580  }
1581  }
1582 
1583 
1584  // edgeProcAddressing
1585  // ~~~~~~~~~~~~~~~~~~
1586  // 1. straight mapping from proc faMesh back to reconstructed
1587  // faMesh
1588  determineFaEdgeMapping
1589  (
1590  masterPatch, // reconstructed faMesh patch
1591  procFaMeshes, // individual faMeshes
1592  pointProcAddressing, // procPolyMesh to reconstructed
1593 
1594  faEdgeProcAddressing
1595  );
1596 
1597 
1598  // Per patch all edges on the masterPatch
1599  labelListList singlePatchEdgeLabels;
1600 
1601  // 2. edgeProcAddressing : fix ordering on patches
1602  sortFaEdgeMapping
1603  (
1604  masterPatch, // reconstructed faMesh patch
1605  procFaMeshes, // individual faMeshes
1606  pointProcAddressing, // procPolyMesh to reconstructed
1607 
1608  faEdgeProcAddressing, // per proc the map to the master
1609  singlePatchEdgeLabels // per patch the map to the master
1610  );
1611 
1612 
1613  const faMesh& procMesh0 = procFaMeshes[0];
1614 
1615 
1616  // Add faPatches
1617  // ~~~~~~~~~~~~~
1618 
1619  faPatchList completePatches(singlePatchEdgeLabels.size());
1620  label nPatches = 0;
1621  forAll(completePatches, patchi)
1622  {
1623  const labelList& patchEdgeLabels =
1624  singlePatchEdgeLabels[patchi];
1625 
1626  const faPatch& fap = procMesh0.boundary()[patchi];
1627 
1628  if (isA<ignoreFaPatch>(fap))
1629  {
1630  // These are not real edges
1631  continue;
1632  }
1633 
1634  const label neiPolyPatchId = fap.ngbPolyPatchIndex();
1635 
1636  completePatches.set
1637  (
1638  nPatches,
1639  fap.clone
1640  (
1641  masterFaMesh.boundary(),
1642  patchEdgeLabels,
1643  nPatches, // index
1644  neiPolyPatchId
1645  )
1646  );
1647 
1648  ++nPatches;
1649  }
1650 
1651  completePatches.resize(nPatches);
1652 
1653  // Serial mesh - no parallel communication
1654 
1655  const bool oldParRun = Pstream::parRun(false);
1656 
1657  masterFaMesh.addFaPatches(completePatches);
1658 
1659  Pstream::parRun(oldParRun); // Restore parallel state
1660 
1661 
1662  // Write mesh & individual addressing
1663  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1664 
1666  (
1667  masterMesh.time().timeName(),
1668  masterFaMesh,
1669  masterFaMesh.faceLabels()
1670  );
1671 
1672  forAll(procFaMeshes, proci)
1673  {
1674  const faMesh& procMesh = procFaMeshes[proci];
1675 
1677  (
1678  IOobject
1679  (
1680  "procAddressing",
1681  masterMesh.time().timeName(),
1683  procMesh.thisDb(),
1687  ),
1688  faBoundProcAddressing[proci],
1689  faFaceProcAddressing[proci],
1690  faPointProcAddressing[proci],
1691  faEdgeProcAddressing[proci]
1692  );
1693  }
1694  }
1695  }
1696  }
1697 
1698  Info<< "End\n" << endl;
1699 
1700  return 0;
1701 }
1702 
1703 
1704 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:202
label nPatches
Definition: readKivaGrid.H:394
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:133
const polyBoundaryMesh & pbm
const Type & value() const noexcept
Return const reference to value.
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
Definition: fvMeshTools.C:317
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:476
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
void writeMesh() const
Write equivalent mesh information at the polyMesh faceInstances time.
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Definition: createK.H:12
A class for handling file names.
Definition: fileName.H:72
readOption readOpt() const noexcept
Get the read option.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:844
label nPoints() const noexcept
Number of mesh points.
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:600
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:526
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:832
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
const fileName & caseName() const noexcept
Return case name.
Definition: TimePathsI.H:78
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
"ascii" (normal default)
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:697
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
engineTime & runTime
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch ("procBoundary..") constructed from the pair of processor IDs...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1586
virtual autoPtr< faPatch > clone(const faBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: faPatch.H:257
fileName relativePath(const fileName &input, const bool caseTag=false) const
Return the input relative to the globalPath by stripping off a leading value of the globalPath...
Definition: TimePathsI.H:103
label nInternalEdges() const
Number of internal edges.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:418
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:388
static void noParallel()
Remove the parallel options.
Definition: argList.C:598
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:168
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
wordList regionNames
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
static void writeContents(const IOobject &io, const UList< T > &content)
Write contents. The IOobject is never registered.
Definition: IOList.C:141
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: refPtrI.H:216
label nFaces() const noexcept
Number of mesh faces.
const labelList & edgeLabels() const noexcept
Return the list of edges.
Definition: faPatch.H:335
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
instantList select(const instantList &times) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:88
bool processorCase() const noexcept
True if this is a processor case.
Definition: TimePathsI.H:52
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2252
const labelList & meshPoints() const
Return labelList of mesh points in patch.
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition: faMesh.H:1203
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
word timeName
Definition: getTimeIndex.H:3
A list of faces which address into the list of points.
bool allowLibs() const
The controlDict &#39;libs&#39; entry is allowed to be used. (eg, has not been disabled by the -no-libs option...
Definition: argList.C:2270
void add(const boundBox &bb)
Extend to include the second box.
Definition: boundBoxI.H:323
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition: topoSet.C:693
const Map< label > & meshPointMap() const
Mesh point map.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:325
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:960
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
word findInstance(const fileName &directory, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null, const bool constant_fallback=true) const
Return time instance (location) of directory containing the file name (eg, used in reading mesh data)...
Definition: Time.C:725
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
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
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: faMesh.C:1015
const Time & time() const noexcept
Return time registry.
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:295
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition: boundBoxI.H:198
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true, const bool fullyMapped=false)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
Definition: fvMeshAdder.C:67
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1101
Reading is optional [identical to LAZY_READ].
static const word null
An empty word.
Definition: word.H:84
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:399
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
label nInternalFaces() const noexcept
Number of internal faces.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:609
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
const fileName & caseName() const noexcept
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:62
static label procPatchPairs(const UPtrList< polyMesh > &meshes, List< DynamicList< label >> &localPatch, List< DynamicList< label >> &remoteMesh, List< DynamicList< label >> &remotePatch)
Helper: find pairs of processor patches. Return number of non-processor patches.
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1068
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:534
label nEdges() const
Number of edges in patch.
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
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
const fileName & rootPath() const noexcept
Return root path.
Definition: argListI.H:56
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:72
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:24
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:74
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeAddressing() const
Write proc addressing at the polyMesh faceInstances time.
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:663
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:750
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:42
static void patchFacePairs(const UPtrList< polyMesh > &meshes, const List< DynamicList< label >> &localPatch, const List< DynamicList< label >> &remoteMesh, const List< DynamicList< label >> &remotePatch, labelListList &localBoundaryFace, labelListList &remoteFaceMesh, labelListList &remoteBoundaryFace)
Helper: expand list of coupled patches into pairs of coupled faces.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
label newPointi
Definition: readKivaGrid.H:497
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
const word & regionDir
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:54
label nNonProcessor() const
The number of patches before the first processor patch.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:826
Required Classes.
Nothing to be read.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
IOstreamOption::streamFormat writeFormat() const noexcept
Get write stream format.
Definition: TimeI.H:123
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
A primitive field of type <T> with automated input and output.
static void writeMaps(Ostream &os, const word &key, const labelListList &maps)
static const word & regionName(const word &region)
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:796
Do not request registration (bool: false)
const dimensionedScalar mp
Proton mass.
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition: HashTableI.H:174
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
label ngbPolyPatchIndex() const noexcept
The neighbour polyPatch index.
Definition: faPatch.H:359
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127