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-2024 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 using namespace Foam;
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
81 // usually meshes get written with limited precision (6 digits)
82 static const scalar defaultMergeTol = 1e-7;
83 
84 
85 // Determine which faces are coupled. Uses geometric merge distance.
86 // Looks either at all boundaryFaces (fullMatch) or only at the
87 // procBoundaries for proci. Assumes that masterMesh contains already merged
88 // all the processors < proci.
89 autoPtr<faceCoupleInfo> determineCoupledFaces
90 (
91  const bool fullMatch,
92  const label masterMeshProcStart,
93  const label masterMeshProcEnd,
94  const polyMesh& masterMesh,
95  const label meshToAddProcStart,
96  const label meshToAddProcEnd,
97  const polyMesh& meshToAdd,
98  const scalar mergeDist
99 )
100 {
101  if (fullMatch || masterMesh.nCells() == 0)
102  {
104  (
105  masterMesh,
106  meshToAdd,
107  mergeDist, // Absolute merging distance
108  true // Matching faces identical
109  );
110  }
111  else
112  {
113  // Pick up all patches on masterMesh ending in "toDDD" where DDD is
114  // the processor number proci.
115 
116  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
117 
118 
119  DynamicList<label> masterFaces
120  (
121  masterMesh.nFaces()
122  - masterMesh.nInternalFaces()
123  );
124 
125 
126  forAll(masterPatches, patchi)
127  {
128  const polyPatch& pp = masterPatches[patchi];
129 
130  if (isA<processorPolyPatch>(pp))
131  {
132  for
133  (
134  label proci=meshToAddProcStart;
135  proci<meshToAddProcEnd;
136  proci++
137  )
138  {
139  const string toProcString("to" + name(proci));
140  if (
141  pp.name().rfind(toProcString)
142  == (pp.name().size()-toProcString.size())
143  )
144  {
145  label meshFacei = pp.start();
146  forAll(pp, i)
147  {
148  masterFaces.append(meshFacei++);
149  }
150  break;
151  }
152  }
153 
154  }
155  }
156  masterFaces.shrink();
157 
158 
159  // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
160  // where DDD is the processor number proci and YYY is < proci.
161 
162  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
163 
164  DynamicList<label> addFaces
165  (
166  meshToAdd.nFaces()
167  - meshToAdd.nInternalFaces()
168  );
169 
170  forAll(addPatches, patchi)
171  {
172  const polyPatch& pp = addPatches[patchi];
173 
174  if (isA<processorPolyPatch>(pp))
175  {
176  bool isConnected = false;
177 
178  for
179  (
180  label mergedProci=masterMeshProcStart;
181  !isConnected && (mergedProci < masterMeshProcEnd);
182  mergedProci++
183  )
184  {
185  for
186  (
187  label proci = meshToAddProcStart;
188  proci < meshToAddProcEnd;
189  proci++
190  )
191  {
192  const word fromProcString
193  (
194  processorPolyPatch::newName(proci, mergedProci)
195  );
196 
197  if (pp.name() == fromProcString)
198  {
199  isConnected = true;
200  break;
201  }
202  }
203  }
204 
205  if (isConnected)
206  {
207  label meshFacei = pp.start();
208  forAll(pp, i)
209  {
210  addFaces.append(meshFacei++);
211  }
212  }
213  }
214  }
215  addFaces.shrink();
216 
218  (
219  masterMesh,
220  masterFaces,
221  meshToAdd,
222  addFaces,
223  mergeDist, // Absolute merging distance
224  true, // Matching faces identical?
225  false, // If perfect match are faces already ordered
226  // (e.g. processor patches)
227  false // are faces each on separate patch?
228  );
229  }
230 }
231 
232 
233 autoPtr<mapPolyMesh> mergeSharedPoints
234 (
235  const scalar mergeDist,
236  polyMesh& mesh,
237  labelListList& pointProcAddressing
238 )
239 {
240  // Find out which sets of points get merged and create a map from
241  // mesh point to unique point.
242  Map<label> pointToMaster
243  (
245  (
246  mesh,
247  mergeDist
248  )
249  );
250 
251  Info<< "mergeSharedPoints : detected " << pointToMaster.size()
252  << " points that are to be merged." << endl;
253 
254  if (returnReduceAnd(pointToMaster.empty()))
255  {
256  return nullptr;
257  }
258 
259  polyTopoChange meshMod(mesh);
260 
261  fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
262 
263  // Change the mesh (no inflation). Note: parallel comms allowed.
264  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
265 
266  // Update fields. No inflation, parallel sync.
267  mesh.updateMesh(map());
268 
269  // pointProcAddressing give indices into the master mesh so adapt them
270  // for changed point numbering.
271 
272  // Adapt constructMaps for merged points.
273  forAll(pointProcAddressing, proci)
274  {
275  labelList& constructMap = pointProcAddressing[proci];
276 
277  forAll(constructMap, i)
278  {
279  label oldPointi = constructMap[i];
280 
281  // New label of point after changeMesh.
282  label newPointi = map().reversePointMap()[oldPointi];
283 
284  if (newPointi < -1)
285  {
286  constructMap[i] = -newPointi-2;
287  }
288  else if (newPointi >= 0)
289  {
290  constructMap[i] = newPointi;
291  }
292  else
293  {
295  << "Problem. oldPointi:" << oldPointi
296  << " newPointi:" << newPointi << abort(FatalError);
297  }
298  }
299  }
300 
301  return map;
302 }
303 
304 
305 boundBox procBounds
306 (
307  const PtrList<Time>& databases,
308  const word& regionName
309 )
310 {
311  boundBox bb;
312 
313  for (const Time& procDb : databases)
314  {
315  fileName pointsInstance
316  (
317  procDb.findInstance(polyMesh::meshDir(regionName), "points")
318  );
319 
320  if (pointsInstance != procDb.timeName())
321  {
323  << "Your time was specified as " << procDb.timeName()
324  << " but there is no polyMesh/points in that time." << nl
325  << "(points file at " << pointsInstance << ')' << nl
326  << "Please rerun with the correct time specified"
327  << " (through the -constant, -time or -latestTime "
328  << "(at your option)."
329  << endl << exit(FatalError);
330  }
331 
332  Info<< "Reading points from "
333  << procDb.caseName()
334  << " for time = " << procDb.timeName()
335  << nl << endl;
336 
338  (
339  IOobject
340  (
341  "points",
342  pointsInstance,
344  procDb,
348  )
349  );
350 
351  bb.add(points);
352  }
353 
354  return bb;
355 }
356 
357 
358 void writeDistribution
359 (
360  Time& runTime,
361  const fvMesh& masterMesh,
362  const labelListList& cellProcAddressing
363 )
364 {
365  // Write the decomposition as labelList for use with 'manual'
366  // decomposition method.
367  labelIOList cellDecomposition
368  (
369  IOobject
370  (
371  "cellDecomposition",
372  masterMesh.facesInstance(),
373  masterMesh,
377  ),
378  masterMesh.nCells()
379  );
380 
381  forAll(cellProcAddressing, proci)
382  {
383  const labelList& pCells = cellProcAddressing[proci];
384  labelUIndList(cellDecomposition, pCells) = proci;
385  }
386 
387  cellDecomposition.write();
388 
389  Info<< nl << "Wrote decomposition to "
390  << cellDecomposition.objectRelPath()
391  << " for use in manual decomposition." << endl;
392 
393  // Write as volScalarField for postprocessing.
394  // Change time to 0 if was 'constant'
395  {
396  const scalar oldTime = runTime.value();
397  const label oldIndex = runTime.timeIndex();
398  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
399  {
400  runTime.setTime(0, oldIndex+1);
401  }
402 
403  volScalarField cellDist
404  (
405  IOobject
406  (
407  "cellDist",
408  runTime.timeName(),
409  masterMesh,
413  ),
414  masterMesh,
417  );
418 
419  forAll(cellDecomposition, celli)
420  {
421  cellDist[celli] = cellDecomposition[celli];
422  }
423 
424  cellDist.correctBoundaryConditions();
425  cellDist.write();
426 
427  Info<< nl << "Wrote decomposition to "
428  << cellDist.objectRelPath()
429  << " (volScalarField) for visualization."
430  << endl;
431 
432  // Restore time
433  runTime.setTime(oldTime, oldIndex);
434  }
435 }
436 
437 
438 void writeMesh
439 (
440  const fvMesh& mesh,
441  const labelListList& cellProcAddressing
442 )
443 {
444  const fileName outputDir
445  (
446  mesh.time().path()
447  / mesh.time().timeName()
448  / mesh.regionName()
450  );
451 
452  Info<< nl << "Writing merged mesh to "
453  << mesh.time().relativePath(outputDir) << nl << endl;
454 
455  if (!mesh.write())
456  {
458  << "Failed writing polyMesh."
459  << exit(FatalError);
460  }
462 }
463 
464 
465 void writeMaps
466 (
467  const label masterInternalFaces,
468  const labelUList& masterOwner,
469  const polyMesh& procMesh,
470  const labelUList& cellProcAddressing,
471  const labelUList& faceProcAddressing,
472  const labelUList& pointProcAddressing,
473  const labelUList& boundProcAddressing
474 )
475 {
476  const fileName outputDir
477  (
478  procMesh.time().caseName()
479  / procMesh.facesInstance()
480  / procMesh.regionName()
482  );
483 
484  IOobject ioAddr
485  (
486  "procAddressing",
487  procMesh.facesInstance(),
489  procMesh,
493  );
494 
495 
496  Info<< "Writing addressing : " << outputDir << nl;
497 
498  // From processor point to reconstructed mesh point
499 
500  Info<< " pointProcAddressing" << endl;
501  ioAddr.rename("pointProcAddressing");
502  IOListRef<label>(ioAddr, pointProcAddressing).write();
503 
504  // From processor face to reconstructed mesh face
505  Info<< " faceProcAddressing" << endl;
506  ioAddr.rename("faceProcAddressing");
507  labelIOList faceProcAddr(ioAddr, faceProcAddressing);
508 
509  // Now add turning index to faceProcAddressing.
510  // See reconstructPar for meaning of turning index.
511  forAll(faceProcAddr, procFacei)
512  {
513  const label masterFacei = faceProcAddr[procFacei];
514 
515  if
516  (
517  !procMesh.isInternalFace(procFacei)
518  && masterFacei < masterInternalFaces
519  )
520  {
521  // proc face is now external but used to be internal face.
522  // Check if we have owner or neighbour.
523 
524  label procOwn = procMesh.faceOwner()[procFacei];
525  label masterOwn = masterOwner[masterFacei];
526 
527  if (cellProcAddressing[procOwn] == masterOwn)
528  {
529  // No turning. Offset by 1.
530  faceProcAddr[procFacei]++;
531  }
532  else
533  {
534  // Turned face.
535  faceProcAddr[procFacei] = -1 - faceProcAddr[procFacei];
536  }
537  }
538  else
539  {
540  // No turning. Offset by 1.
541  faceProcAddr[procFacei]++;
542  }
543  }
544 
545  faceProcAddr.write();
546 
547 
548  // From processor cell to reconstructed mesh cell
549  Info<< " cellProcAddressing" << endl;
550  ioAddr.rename("cellProcAddressing");
551  IOListRef<label>(ioAddr, cellProcAddressing).write();
552 
553 
554  // From processor patch to reconstructed mesh patch
555  Info<< " boundaryProcAddressing" << endl;
556  ioAddr.rename("boundaryProcAddressing");
557  IOListRef<label>(ioAddr, boundProcAddressing).write();
558 
559  Info<< endl;
560 }
561 
562 
563 void printWarning()
564 {
565  Info<<
566 "Merge individual processor meshes back into one master mesh.\n"
567 "Use if the original master mesh has been deleted or the processor meshes\n"
568 "have been modified (topology change).\n"
569 "This tool will write the resulting mesh to a new time step and construct\n"
570 "xxxxProcAddressing files in the processor meshes so reconstructPar can be\n"
571 "used to regenerate the fields on the master mesh.\n\n"
572 "Not well tested & use at your own risk!\n\n";
573 }
574 
575 
576 int main(int argc, char *argv[])
577 {
579  (
580  "Reconstruct a mesh using geometric/topological information only"
581  );
582 
583  // Enable -constant ... if someone really wants it
584  // Enable -withZero to prevent accidentally trashing the initial fields
585  timeSelector::addOptions(true, true); // constant(true), zero(true)
586 
588 
591  (
592  "addressing-only",
593  "Create procAddressing only without overwriting the mesh"
594  );
596  (
597  "mergeTol",
598  "scalar",
599  "The merge distance relative to the bounding box size (default 1e-7)"
600  );
602  (
603  "fullMatch",
604  "Do (slower) geometric matching on all boundary faces"
605  );
607  (
608  "procMatch",
609  "Do matching on processor faces only"
610  );
612  (
613  "cellDist",
614  "Write cell distribution as a labelList - for use with 'manual' "
615  "decomposition method or as a volScalarField for post-processing."
616  );
617 
618  #include "addAllRegionOptions.H"
619 
620  #include "setRootCase.H"
621  #include "createTime.H"
622 
623  printWarning();
624 
625  const bool fullMatch = args.found("fullMatch");
626  const bool procMatch = args.found("procMatch");
627  const bool writeCellDist = args.found("cellDist");
628 
629  const bool writeAddrOnly = args.found("addressing-only");
630 
631  const scalar mergeTol =
632  args.getOrDefault<scalar>("mergeTol", defaultMergeTol);
633 
634  if (fullMatch)
635  {
636  Info<< "Use geometric matching on all boundary faces." << nl << endl;
637  }
638  else if (procMatch)
639  {
640  Info<< "Use geometric matching on correct procBoundaries only." << nl
641  << "This assumes a correct decomposition." << endl;
642  }
643  else
644  {
645  Info<< "Merge assuming correct, fully matched procBoundaries." << nl
646  << endl;
647  }
648 
649  if (fullMatch || procMatch)
650  {
651  const scalar writeTol =
652  std::pow(scalar(10), -scalar(IOstream::defaultPrecision()));
653 
654  Info<< "Merge tolerance : " << mergeTol << nl
655  << "Write tolerance : " << writeTol << endl;
656 
657  if
658  (
660  && mergeTol < writeTol
661  )
662  {
664  << "Your current settings specify ASCII writing with "
665  << IOstream::defaultPrecision() << " digits precision." << endl
666  << "Your merging tolerance (" << mergeTol << ")"
667  << " is finer than this." << endl
668  << "Please change your writeFormat to binary"
669  << " or increase the writePrecision" << endl
670  << "or adjust the merge tolerance (-mergeTol)."
671  << exit(FatalError);
672  }
673  }
674 
675  // Get region names
676  #include "getAllRegionOptions.H"
677 
678  // Determine the processor count
679  label nProcs{0};
680 
681  if (regionNames.empty())
682  {
684  << "No regions specified or detected."
685  << exit(FatalError);
686  }
687  else if (regionNames[0] == polyMesh::defaultRegion)
688  {
689  nProcs = fileHandler().nProcs(args.path());
690  }
691  else
692  {
693  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
694 
695  if (regionNames.size() == 1)
696  {
697  Info<< "Using region: " << regionNames[0] << nl << endl;
698  }
699  }
700 
701  if (!nProcs)
702  {
704  << "No processor* directories found"
705  << exit(FatalError);
706  }
707 
708  // Read all time databases
709  PtrList<Time> databases(nProcs);
710 
711  Info<< "Found " << nProcs << " processor directories" << endl;
712  forAll(databases, proci)
713  {
714  Info<< " Reading database "
715  << args.caseName()/("processor" + Foam::name(proci))
716  << endl;
717 
718  databases.set
719  (
720  proci,
721  new Time
722  (
724  args.rootPath(),
725  args.caseName()/("processor" + Foam::name(proci)),
727  args.allowLibs()
728  )
729  );
730  }
731  Info<< endl;
732 
733  // Use the times list from the master processor
734  // and select a subset based on the command-line options
736  (
737  databases[0].times(),
738  args
739  );
740 
741  // Loop over all times
742  forAll(timeDirs, timei)
743  {
744  // Set time for global database
745  runTime.setTime(timeDirs[timei], timei);
746 
747  // Set time for all databases
748  forAll(databases, proci)
749  {
750  databases[proci].setTime(timeDirs[timei], timei);
751  }
752 
753  Info<< "Time = " << runTime.timeName() << endl;
754 
755  // Check for any mesh changes
756  label nMeshChanged = 0;
757  boolList hasRegionMesh(regionNames.size(), false);
758  forAll(regionNames, regioni)
759  {
760  const word& regionName = regionNames[regioni];
761 
762  IOobject facesIO
763  (
764  "faces",
765  databases[0].timeName(),
767  databases[0],
770  );
771 
772  // Problem: faceCompactIOList recognises both 'faceList' and
773  // 'faceCompactList' so we should be lenient when doing
774  // typeHeaderOk
775 
776  const bool ok = facesIO.typeHeaderOk<faceCompactIOList>(false);
777 
778  if (ok)
779  {
780  hasRegionMesh[regioni] = true;
781  ++nMeshChanged;
782  }
783  }
784 
785  // Check for any mesh changes
786  if (!nMeshChanged)
787  {
788  Info<< "No meshes" << nl << endl;
789  continue;
790  }
791 
792  Info<< endl;
793 
794  forAll(regionNames, regioni)
795  {
796  const word& regionName = regionNames[regioni];
798 
799  if (!hasRegionMesh[regioni])
800  {
801  Info<< "region=" << regionName << " (no mesh)" << nl << endl;
802  continue;
803  }
804 
805  if (regionNames.size() > 1)
806  {
807  Info<< "region=" << regionName << nl;
808  }
809 
810 
811  // Addressing from processor to reconstructed case
812  labelListList cellProcAddressing(nProcs);
813  labelListList faceProcAddressing(nProcs);
814  labelListList pointProcAddressing(nProcs);
815  labelListList boundProcAddressing(nProcs);
816 
817 
818  // Internal faces on the final reconstructed mesh
819  label masterInternalFaces;
820 
821  // Owner addressing on the final reconstructed mesh
822  labelList masterOwner;
823 
824  if (procMatch)
825  {
826  // Read point on individual processors to determine
827  // merge tolerance
828  // (otherwise single cell domains might give problems)
829 
830  const boundBox bb = procBounds(databases, regionDir);
831  const scalar mergeDist = mergeTol*bb.mag();
832 
833  Info<< "Overall mesh bounding box : " << bb << nl
834  << "Relative tolerance : " << mergeTol << nl
835  << "Absolute matching distance : " << mergeDist << nl
836  << endl;
837 
838 
839  // Construct empty mesh.
840  PtrList<fvMesh> masterMesh(nProcs);
841 
842  for (label proci=0; proci<nProcs; proci++)
843  {
844  masterMesh.set
845  (
846  proci,
847  new fvMesh
848  (
849  IOobject
850  (
851  regionName,
852  runTime.timeName(),
853  runTime,
855  ),
856  Zero
857  )
858  );
859 
860  fvMesh meshToAdd
861  (
862  IOobject
863  (
864  regionName,
865  databases[proci].timeName(),
866  databases[proci]
867  )
868  );
869 
870  // Initialize its addressing
871  cellProcAddressing[proci] = identity(meshToAdd.nCells());
872  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
873  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
874  boundProcAddressing[proci] =
875  identity(meshToAdd.boundaryMesh().size());
876 
877  // Find geometrically shared points/faces.
878  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
879  (
880  fullMatch,
881  proci,
882  proci,
883  masterMesh[proci],
884  proci,
885  proci,
886  meshToAdd,
887  mergeDist
888  );
889 
890  // Add elements to mesh
892  (
893  masterMesh[proci],
894  meshToAdd,
895  couples()
896  );
897 
898  // Added processor
899  renumber(map().addedCellMap(), cellProcAddressing[proci]);
900  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
901  renumber(map().addedPointMap(), pointProcAddressing[proci]);
902  renumber(map().addedPatchMap(), boundProcAddressing[proci]);
903  }
904  for (label step=2; step<nProcs*2; step*=2)
905  {
906  for (label proci=0; proci<nProcs; proci+=step)
907  {
908  label next = proci + step/2;
909  if(next >= nProcs)
910  {
911  continue;
912  }
913 
914  Info<< "Merging mesh " << proci << " with "
915  << next << endl;
916 
917  // Find geometrically shared points/faces.
918  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
919  (
920  fullMatch,
921  proci,
922  next,
923  masterMesh[proci],
924  next,
925  proci+step,
926  masterMesh[next],
927  mergeDist
928  );
929 
930  // Add elements to mesh
932  (
933  masterMesh[proci],
934  masterMesh[next],
935  couples()
936  );
937 
938  // Processors that were already in masterMesh
939  for (label mergedI=proci; mergedI<next; mergedI++)
940  {
941  renumber
942  (
943  map().oldCellMap(),
944  cellProcAddressing[mergedI]
945  );
946 
947  renumber
948  (
949  map().oldFaceMap(),
950  faceProcAddressing[mergedI]
951  );
952 
953  renumber
954  (
955  map().oldPointMap(),
956  pointProcAddressing[mergedI]
957  );
958 
959  // Note: boundary is special since can contain -1.
960  renumber
961  (
962  map().oldPatchMap(),
963  boundProcAddressing[mergedI]
964  );
965  }
966 
967  // Added processor
968  for
969  (
970  label addedI=next;
971  addedI<min(proci+step, nProcs);
972  addedI++
973  )
974  {
975  renumber
976  (
977  map().addedCellMap(),
978  cellProcAddressing[addedI]
979  );
980 
981  renumber
982  (
983  map().addedFaceMap(),
984  faceProcAddressing[addedI]
985  );
986 
987  renumber
988  (
989  map().addedPointMap(),
990  pointProcAddressing[addedI]
991  );
992 
993  renumber
994  (
995  map().addedPatchMap(),
996  boundProcAddressing[addedI]
997  );
998  }
999 
1000  masterMesh.set(next, nullptr);
1001  }
1002  }
1003 
1004  for (label proci=0; proci<nProcs; proci++)
1005  {
1006  Info<< "Reading mesh to add from "
1007  << databases[proci].caseName()
1008  << " for time = " << databases[proci].timeName()
1009  << nl << nl << endl;
1010  }
1011 
1012  // See if any points on the mastermesh have become connected
1013  // because of connections through processor meshes.
1014  mergeSharedPoints(mergeDist,masterMesh[0],pointProcAddressing);
1015 
1016  // Save some properties on the reconstructed mesh
1017  masterInternalFaces = masterMesh[0].nInternalFaces();
1018  masterOwner = masterMesh[0].faceOwner();
1019 
1020 
1021  if (writeAddrOnly)
1022  {
1023  Info<< nl
1024  << "Disabled writing of merged mesh (-addressing-only)"
1025  << nl << nl;
1026  }
1027  else
1028  {
1029  // Write reconstructed mesh
1030  writeMesh(masterMesh[0], cellProcAddressing);
1031  if (writeCellDist)
1032  {
1033  writeDistribution
1034  (
1035  runTime,
1036  masterMesh[0],
1037  cellProcAddressing
1038  );
1039  }
1040  }
1041  }
1042  else
1043  {
1044  // Load all meshes
1045  PtrList<fvMesh> fvMeshes(nProcs);
1046  for (label proci=0; proci<nProcs; proci++)
1047  {
1048  fvMeshes.set
1049  (
1050  proci,
1051  new fvMesh
1052  (
1053  IOobject
1054  (
1055  regionName,
1056  databases[proci].timeName(),
1057  databases[proci]
1058  )
1059  )
1060  );
1061  }
1062 
1063  // Construct pointers to meshes
1064  UPtrList<polyMesh> meshes(fvMeshes.size());
1065  forAll(fvMeshes, proci)
1066  {
1067  meshes.set(proci, &fvMeshes[proci]);
1068  }
1069 
1070  // Get pairs of patches to stitch. These pairs have to
1071  // - have ordered, opposite faces (so one to one correspondence)
1072  List<DynamicList<label>> localPatch;
1073  List<DynamicList<label>> remoteProc;
1074  List<DynamicList<label>> remotePatch;
1075  const label nGlobalPatches = polyMeshAdder::procPatchPairs
1076  (
1077  meshes,
1078  localPatch,
1079  remoteProc,
1080  remotePatch
1081  );
1082 
1083  // Collect matching boundary faces on patches-to-stitch
1084  labelListList localBoundaryFace;
1085  labelListList remoteFaceProc;
1086  labelListList remoteBoundaryFace;
1088  (
1089  meshes,
1090  localPatch,
1091  remoteProc,
1092  remotePatch,
1093  localBoundaryFace,
1094  remoteFaceProc,
1095  remoteBoundaryFace
1096  );
1097 
1098  // All matched faces assumed to have vertex0 matched
1099  labelListList remoteFaceStart(meshes.size());
1100  forAll(meshes, proci)
1101  {
1102  const labelList& procFaces = localBoundaryFace[proci];
1103  remoteFaceStart[proci].setSize(procFaces.size(), 0);
1104  }
1105 
1106 
1107 
1108  labelListList patchMap(meshes.size());
1109  labelListList pointZoneMap(meshes.size());
1110  labelListList faceZoneMap(meshes.size());
1111  labelListList cellZoneMap(meshes.size());
1112  forAll(meshes, proci)
1113  {
1114  const polyMesh& mesh = meshes[proci];
1115  patchMap[proci] = identity(mesh.boundaryMesh().size());
1116 
1117  // Remove excess patches
1118  patchMap[proci].setSize(nGlobalPatches);
1119 
1120  pointZoneMap[proci] = identity(mesh.pointZones().size());
1121  faceZoneMap[proci] = identity(mesh.faceZones().size());
1122  cellZoneMap[proci] = identity(mesh.cellZones().size());
1123  }
1124 
1125 
1126  refPtr<fvMesh> masterMeshPtr;
1127  {
1128  // Do in-place addition on proc0.
1129 
1130  const labelList oldFaceOwner(fvMeshes[0].faceOwner());
1131 
1133  (
1134  0, // index of mesh to modify (== mesh_)
1135  fvMeshes,
1136  oldFaceOwner,
1137 
1138  // Coupling info
1139  localBoundaryFace,
1140  remoteFaceProc,
1141  remoteBoundaryFace,
1142 
1143  boundProcAddressing,
1144  cellProcAddressing,
1145  faceProcAddressing,
1146  pointProcAddressing
1147  );
1148 
1149  // Remove zero-faces processor patches
1150  const polyBoundaryMesh& pbm = fvMeshes[0].boundaryMesh();
1151  labelList oldToNew(pbm.size(), -1);
1152  label newi = 0;
1153  // Non processor patches first
1154  forAll(pbm, patchi)
1155  {
1156  const auto& pp = pbm[patchi];
1157  if (!isA<processorPolyPatch>(pp) || pp.size())
1158  {
1159  oldToNew[patchi] = newi++;
1160  }
1161  }
1162  const label nNonProcPatches = newi;
1163 
1164  // Move all deletable patches to the end
1165  forAll(oldToNew, patchi)
1166  {
1167  if (oldToNew[patchi] == -1)
1168  {
1169  oldToNew[patchi] = newi++;
1170  }
1171  }
1173  (
1174  fvMeshes[0],
1175  oldToNew,
1176  nNonProcPatches,
1177  false
1178  );
1179 
1180  masterMeshPtr.cref(fvMeshes[0]);
1181  }
1182 
1183 
1184  const fvMesh& masterMesh = masterMeshPtr();
1185 
1186  // Number of internal faces on the final reconstructed mesh
1187  masterInternalFaces = masterMesh.nInternalFaces();
1188 
1189  // Owner addressing on the final reconstructed mesh
1190  masterOwner = masterMesh.faceOwner();
1191 
1192 
1193  // Write reconstructed mesh
1194  // Override:
1195  // - caseName
1196  // - processorCase flag
1197  // so the resulting mesh goes to the correct location (even with
1198  // collated). The better way of solving this is to construct
1199  // (zero) mesh on the undecomposed runTime.
1200 
1201  if (writeAddrOnly)
1202  {
1203  Info<< nl
1204  << "Disabled writing of merged mesh (-addressing-only)"
1205  << nl << nl;
1206  }
1207  else
1208  {
1209  Time& masterTime = const_cast<Time&>(masterMesh.time());
1210 
1211  const word oldCaseName = masterTime.caseName();
1212  masterTime.caseName() = runTime.caseName();
1213  const bool oldProcCase(masterTime.processorCase(false));
1214 
1215  // Write reconstructed mesh
1216  writeMesh(masterMesh, cellProcAddressing);
1217  if (writeCellDist)
1218  {
1219  writeDistribution
1220  (
1221  runTime,
1222  masterMesh,
1223  cellProcAddressing
1224  );
1225  }
1226  masterTime.caseName() = oldCaseName;
1227  masterTime.processorCase(oldProcCase);
1228  }
1229  }
1230 
1231 
1232  // Write the addressing
1233 
1234  Info<< "Reconstructing addressing from processor meshes"
1235  << " to the newly reconstructed mesh" << nl << endl;
1236 
1237  forAll(databases, proci)
1238  {
1239  Info<< "Processor " << proci << nl
1240  << "Read processor mesh: "
1241  << (databases[proci].caseName()/regionDir) << endl;
1242 
1243  polyMesh procMesh
1244  (
1245  IOobject
1246  (
1247  regionName,
1248  databases[proci].timeName(),
1249  databases[proci]
1250  )
1251  );
1252 
1253  writeMaps
1254  (
1255  masterInternalFaces,
1256  masterOwner,
1257  procMesh,
1258  cellProcAddressing[proci],
1259  faceProcAddressing[proci],
1260  pointProcAddressing[proci],
1261  boundProcAddressing[proci]
1262  );
1263  }
1264  }
1265  }
1266 
1267  Info<< "End\n" << endl;
1268 
1269  return 0;
1270 }
1271 
1272 
1273 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
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:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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
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:859
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:608
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
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:675
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:531
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
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:374
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
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
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.
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:2159
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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)
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:2177
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
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:320
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
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
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
const Time & time() const noexcept
Return time registry.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
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:1116
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:385
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
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:1113
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static void addVerboseOption(const string &usage="", bool advanced=false)
Enable a &#39;verbose&#39; bool option, with usage information.
Definition: argList.C:520
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
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:484
fileName path() const
Return the full path to the (processor local) case.
Definition: argListI.H:74
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:663
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:496
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:55
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:841
Required Classes.
Nothing to be read.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
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
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)
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: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:811
Do not request registration (bool: false)
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.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127