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-2023 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
318  (
320  "points"
321  )
322  );
323 
324  if (pointsInstance != procDb.timeName())
325  {
327  << "Your time was specified as " << procDb.timeName()
328  << " but there is no polyMesh/points in that time." << nl
329  << "(points file at " << pointsInstance << ')' << nl
330  << "Please rerun with the correct time specified"
331  << " (through the -constant, -time or -latestTime "
332  << "(at your option)."
333  << endl << exit(FatalError);
334  }
335 
336  Info<< "Reading points from "
337  << procDb.caseName()
338  << " for time = " << procDb.timeName()
339  << nl << endl;
340 
342  (
343  IOobject
344  (
345  "points",
346  procDb.findInstance
347  (
349  "points"
350  ),
352  procDb,
356  )
357  );
358 
359  bb.add(points);
360  }
361 
362  return bb;
363 }
364 
365 
366 void writeDistribution
367 (
368  Time& runTime,
369  const fvMesh& masterMesh,
370  const labelListList& cellProcAddressing
371 )
372 {
373  // Write the decomposition as labelList for use with 'manual'
374  // decomposition method.
375  labelIOList cellDecomposition
376  (
377  IOobject
378  (
379  "cellDecomposition",
380  masterMesh.facesInstance(),
381  masterMesh,
385  ),
386  masterMesh.nCells()
387  );
388 
389  forAll(cellProcAddressing, proci)
390  {
391  const labelList& pCells = cellProcAddressing[proci];
392  labelUIndList(cellDecomposition, pCells) = proci;
393  }
394 
395  cellDecomposition.write();
396 
397  Info<< nl << "Wrote decomposition to "
398  << cellDecomposition.objectRelPath()
399  << " for use in manual decomposition." << endl;
400 
401  // Write as volScalarField for postprocessing.
402  // Change time to 0 if was 'constant'
403  {
404  const scalar oldTime = runTime.value();
405  const label oldIndex = runTime.timeIndex();
406  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
407  {
408  runTime.setTime(0, oldIndex+1);
409  }
410 
411  volScalarField cellDist
412  (
413  IOobject
414  (
415  "cellDist",
416  runTime.timeName(),
417  masterMesh,
421  ),
422  masterMesh,
425  );
426 
427  forAll(cellDecomposition, celli)
428  {
429  cellDist[celli] = cellDecomposition[celli];
430  }
431 
432  cellDist.correctBoundaryConditions();
433  cellDist.write();
434 
435  Info<< nl << "Wrote decomposition to "
436  << cellDist.objectRelPath()
437  << " (volScalarField) for visualization."
438  << endl;
439 
440  // Restore time
441  runTime.setTime(oldTime, oldIndex);
442  }
443 }
444 
445 
446 void writeMesh
447 (
448  const fvMesh& mesh,
449  const labelListList& cellProcAddressing
450 )
451 {
452  const fileName outputDir
453  (
454  mesh.time().path()
455  / mesh.time().timeName()
456  / mesh.regionName()
458  );
459 
460  Info<< nl << "Writing merged mesh to "
461  << mesh.time().relativePath(outputDir) << nl << endl;
462 
463  if (!mesh.write())
464  {
466  << "Failed writing polyMesh."
467  << exit(FatalError);
468  }
470 }
471 
472 
473 void writeMaps
474 (
475  const label masterInternalFaces,
476  const labelUList& masterOwner,
477  const polyMesh& procMesh,
478  const labelUList& cellProcAddressing,
479  const labelUList& faceProcAddressing,
480  const labelUList& pointProcAddressing,
481  const labelUList& boundProcAddressing
482 )
483 {
484  const fileName outputDir
485  (
486  procMesh.time().caseName()
487  / procMesh.facesInstance()
488  / procMesh.regionName()
490  );
491 
492  IOobject ioAddr
493  (
494  "procAddressing",
495  procMesh.facesInstance(),
497  procMesh,
501  );
502 
503 
504  Info<< "Writing addressing : " << outputDir << nl;
505 
506  // From processor point to reconstructed mesh point
507 
508  Info<< " pointProcAddressing" << endl;
509  ioAddr.rename("pointProcAddressing");
510  IOListRef<label>(ioAddr, pointProcAddressing).write();
511 
512  // From processor face to reconstructed mesh face
513  Info<< " faceProcAddressing" << endl;
514  ioAddr.rename("faceProcAddressing");
515  labelIOList faceProcAddr(ioAddr, faceProcAddressing);
516 
517  // Now add turning index to faceProcAddressing.
518  // See reconstructPar for meaning of turning index.
519  forAll(faceProcAddr, procFacei)
520  {
521  const label masterFacei = faceProcAddr[procFacei];
522 
523  if
524  (
525  !procMesh.isInternalFace(procFacei)
526  && masterFacei < masterInternalFaces
527  )
528  {
529  // proc face is now external but used to be internal face.
530  // Check if we have owner or neighbour.
531 
532  label procOwn = procMesh.faceOwner()[procFacei];
533  label masterOwn = masterOwner[masterFacei];
534 
535  if (cellProcAddressing[procOwn] == masterOwn)
536  {
537  // No turning. Offset by 1.
538  faceProcAddr[procFacei]++;
539  }
540  else
541  {
542  // Turned face.
543  faceProcAddr[procFacei] = -1 - faceProcAddr[procFacei];
544  }
545  }
546  else
547  {
548  // No turning. Offset by 1.
549  faceProcAddr[procFacei]++;
550  }
551  }
552 
553  faceProcAddr.write();
554 
555 
556  // From processor cell to reconstructed mesh cell
557  Info<< " cellProcAddressing" << endl;
558  ioAddr.rename("cellProcAddressing");
559  IOListRef<label>(ioAddr, cellProcAddressing).write();
560 
561 
562  // From processor patch to reconstructed mesh patch
563  Info<< " boundaryProcAddressing" << endl;
564  ioAddr.rename("boundaryProcAddressing");
565  IOListRef<label>(ioAddr, boundProcAddressing).write();
566 
567  Info<< endl;
568 }
569 
570 
571 void printWarning()
572 {
573  Info<<
574 "Merge individual processor meshes back into one master mesh.\n"
575 "Use if the original master mesh has been deleted or the processor meshes\n"
576 "have been modified (topology change).\n"
577 "This tool will write the resulting mesh to a new time step and construct\n"
578 "xxxxProcAddressing files in the processor meshes so reconstructPar can be\n"
579 "used to regenerate the fields on the master mesh.\n\n"
580 "Not well tested & use at your own risk!\n\n";
581 }
582 
583 
584 int main(int argc, char *argv[])
585 {
587  (
588  "Reconstruct a mesh using geometric/topological information only"
589  );
590 
591  // Enable -constant ... if someone really wants it
592  // Enable -withZero to prevent accidentally trashing the initial fields
593  timeSelector::addOptions(true, true); // constant(true), zero(true)
594 
596 
599  (
600  "addressing-only",
601  "Create procAddressing only without overwriting the mesh"
602  );
604  (
605  "mergeTol",
606  "scalar",
607  "The merge distance relative to the bounding box size (default 1e-7)"
608  );
610  (
611  "fullMatch",
612  "Do (slower) geometric matching on all boundary faces"
613  );
615  (
616  "procMatch",
617  "Do matching on processor faces only"
618  );
620  (
621  "cellDist",
622  "Write cell distribution as a labelList - for use with 'manual' "
623  "decomposition method or as a volScalarField for post-processing."
624  );
625 
626  #include "addAllRegionOptions.H"
627 
628  #include "setRootCase.H"
629  #include "createTime.H"
630 
631  printWarning();
632 
633  const bool fullMatch = args.found("fullMatch");
634  const bool procMatch = args.found("procMatch");
635  const bool writeCellDist = args.found("cellDist");
636 
637  const bool writeAddrOnly = args.found("addressing-only");
638 
639  const scalar mergeTol =
640  args.getOrDefault<scalar>("mergeTol", defaultMergeTol);
641 
642  if (fullMatch)
643  {
644  Info<< "Use geometric matching on all boundary faces." << nl << endl;
645  }
646  else if (procMatch)
647  {
648  Info<< "Use geometric matching on correct procBoundaries only." << nl
649  << "This assumes a correct decomposition." << endl;
650  }
651  else
652  {
653  Info<< "Merge assuming correct, fully matched procBoundaries." << nl
654  << endl;
655  }
656 
657  if (fullMatch || procMatch)
658  {
659  const scalar writeTol =
660  Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
661 
662  Info<< "Merge tolerance : " << mergeTol << nl
663  << "Write tolerance : " << writeTol << endl;
664 
665  if
666  (
668  && mergeTol < writeTol
669  )
670  {
672  << "Your current settings specify ASCII writing with "
673  << IOstream::defaultPrecision() << " digits precision." << endl
674  << "Your merging tolerance (" << mergeTol << ")"
675  << " is finer than this." << endl
676  << "Please change your writeFormat to binary"
677  << " or increase the writePrecision" << endl
678  << "or adjust the merge tolerance (-mergeTol)."
679  << exit(FatalError);
680  }
681  }
682 
683  // Get region names
684  #include "getAllRegionOptions.H"
685 
686  // Determine the processor count
687  label nProcs{0};
688 
689  if (regionNames.empty())
690  {
692  << "No regions specified or detected."
693  << exit(FatalError);
694  }
695  else if (regionNames[0] == polyMesh::defaultRegion)
696  {
697  nProcs = fileHandler().nProcs(args.path());
698  }
699  else
700  {
701  nProcs = fileHandler().nProcs(args.path(), regionNames[0]);
702 
703  if (regionNames.size() == 1)
704  {
705  Info<< "Using region: " << regionNames[0] << nl << endl;
706  }
707  }
708 
709  if (!nProcs)
710  {
712  << "No processor* directories found"
713  << exit(FatalError);
714  }
715 
716  // Read all time databases
717  PtrList<Time> databases(nProcs);
718 
719  Info<< "Found " << nProcs << " processor directories" << endl;
720  forAll(databases, proci)
721  {
722  Info<< " Reading database "
723  << args.caseName()/("processor" + Foam::name(proci))
724  << endl;
725 
726  databases.set
727  (
728  proci,
729  new Time
730  (
732  args.rootPath(),
733  args.caseName()/("processor" + Foam::name(proci)),
735  args.allowLibs()
736  )
737  );
738  }
739  Info<< endl;
740 
741  // Use the times list from the master processor
742  // and select a subset based on the command-line options
744  (
745  databases[0].times(),
746  args
747  );
748 
749  // Loop over all times
750  forAll(timeDirs, timei)
751  {
752  // Set time for global database
753  runTime.setTime(timeDirs[timei], timei);
754 
755  // Set time for all databases
756  forAll(databases, proci)
757  {
758  databases[proci].setTime(timeDirs[timei], timei);
759  }
760 
761  Info<< "Time = " << runTime.timeName() << endl;
762 
763  // Check for any mesh changes
764  label nMeshChanged = 0;
765  boolList hasRegionMesh(regionNames.size(), false);
766  forAll(regionNames, regioni)
767  {
768  const word& regionName = regionNames[regioni];
769 
770  IOobject facesIO
771  (
772  "faces",
773  databases[0].timeName(),
775  databases[0],
778  );
779 
780  // Problem: faceCompactIOList recognises both 'faceList' and
781  // 'faceCompactList' so we should be lenient when doing
782  // typeHeaderOk
783 
784  const bool ok = facesIO.typeHeaderOk<faceCompactIOList>(false);
785 
786  if (ok)
787  {
788  hasRegionMesh[regioni] = true;
789  ++nMeshChanged;
790  }
791  }
792 
793  // Check for any mesh changes
794  if (!nMeshChanged)
795  {
796  Info<< "No meshes" << nl << endl;
797  continue;
798  }
799 
800  Info<< endl;
801 
802  forAll(regionNames, regioni)
803  {
804  const word& regionName = regionNames[regioni];
806 
807  if (!hasRegionMesh[regioni])
808  {
809  Info<< "region=" << regionName << " (no mesh)" << nl << endl;
810  continue;
811  }
812 
813  if (regionNames.size() > 1)
814  {
815  Info<< "region=" << regionName << nl;
816  }
817 
818 
819  // Addressing from processor to reconstructed case
820  labelListList cellProcAddressing(nProcs);
821  labelListList faceProcAddressing(nProcs);
822  labelListList pointProcAddressing(nProcs);
823  labelListList boundProcAddressing(nProcs);
824 
825 
826  // Internal faces on the final reconstructed mesh
827  label masterInternalFaces;
828 
829  // Owner addressing on the final reconstructed mesh
830  labelList masterOwner;
831 
832  if (procMatch)
833  {
834  // Read point on individual processors to determine
835  // merge tolerance
836  // (otherwise single cell domains might give problems)
837 
838  const boundBox bb = procBounds(databases, regionDir);
839  const scalar mergeDist = mergeTol*bb.mag();
840 
841  Info<< "Overall mesh bounding box : " << bb << nl
842  << "Relative tolerance : " << mergeTol << nl
843  << "Absolute matching distance : " << mergeDist << nl
844  << endl;
845 
846 
847  // Construct empty mesh.
848  PtrList<fvMesh> masterMesh(nProcs);
849 
850  for (label proci=0; proci<nProcs; proci++)
851  {
852  masterMesh.set
853  (
854  proci,
855  new fvMesh
856  (
857  IOobject
858  (
859  regionName,
860  runTime.timeName(),
861  runTime,
863  ),
864  Zero
865  )
866  );
867 
868  fvMesh meshToAdd
869  (
870  IOobject
871  (
872  regionName,
873  databases[proci].timeName(),
874  databases[proci]
875  )
876  );
877 
878  // Initialize its addressing
879  cellProcAddressing[proci] = identity(meshToAdd.nCells());
880  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
881  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
882  boundProcAddressing[proci] =
883  identity(meshToAdd.boundaryMesh().size());
884 
885  // Find geometrically shared points/faces.
886  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
887  (
888  fullMatch,
889  proci,
890  proci,
891  masterMesh[proci],
892  proci,
893  proci,
894  meshToAdd,
895  mergeDist
896  );
897 
898  // Add elements to mesh
900  (
901  masterMesh[proci],
902  meshToAdd,
903  couples()
904  );
905 
906  // Added processor
907  renumber(map().addedCellMap(), cellProcAddressing[proci]);
908  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
909  renumber(map().addedPointMap(), pointProcAddressing[proci]);
910  renumber(map().addedPatchMap(), boundProcAddressing[proci]);
911  }
912  for (label step=2; step<nProcs*2; step*=2)
913  {
914  for (label proci=0; proci<nProcs; proci+=step)
915  {
916  label next = proci + step/2;
917  if(next >= nProcs)
918  {
919  continue;
920  }
921 
922  Info<< "Merging mesh " << proci << " with "
923  << next << endl;
924 
925  // Find geometrically shared points/faces.
926  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
927  (
928  fullMatch,
929  proci,
930  next,
931  masterMesh[proci],
932  next,
933  proci+step,
934  masterMesh[next],
935  mergeDist
936  );
937 
938  // Add elements to mesh
940  (
941  masterMesh[proci],
942  masterMesh[next],
943  couples()
944  );
945 
946  // Processors that were already in masterMesh
947  for (label mergedI=proci; mergedI<next; mergedI++)
948  {
949  renumber
950  (
951  map().oldCellMap(),
952  cellProcAddressing[mergedI]
953  );
954 
955  renumber
956  (
957  map().oldFaceMap(),
958  faceProcAddressing[mergedI]
959  );
960 
961  renumber
962  (
963  map().oldPointMap(),
964  pointProcAddressing[mergedI]
965  );
966 
967  // Note: boundary is special since can contain -1.
968  renumber
969  (
970  map().oldPatchMap(),
971  boundProcAddressing[mergedI]
972  );
973  }
974 
975  // Added processor
976  for
977  (
978  label addedI=next;
979  addedI<min(proci+step, nProcs);
980  addedI++
981  )
982  {
983  renumber
984  (
985  map().addedCellMap(),
986  cellProcAddressing[addedI]
987  );
988 
989  renumber
990  (
991  map().addedFaceMap(),
992  faceProcAddressing[addedI]
993  );
994 
995  renumber
996  (
997  map().addedPointMap(),
998  pointProcAddressing[addedI]
999  );
1000 
1001  renumber
1002  (
1003  map().addedPatchMap(),
1004  boundProcAddressing[addedI]
1005  );
1006  }
1007 
1008  masterMesh.set(next, nullptr);
1009  }
1010  }
1011 
1012  for (label proci=0; proci<nProcs; proci++)
1013  {
1014  Info<< "Reading mesh to add from "
1015  << databases[proci].caseName()
1016  << " for time = " << databases[proci].timeName()
1017  << nl << nl << endl;
1018  }
1019 
1020  // See if any points on the mastermesh have become connected
1021  // because of connections through processor meshes.
1022  mergeSharedPoints(mergeDist,masterMesh[0],pointProcAddressing);
1023 
1024  // Save some properties on the reconstructed mesh
1025  masterInternalFaces = masterMesh[0].nInternalFaces();
1026  masterOwner = masterMesh[0].faceOwner();
1027 
1028 
1029  if (writeAddrOnly)
1030  {
1031  Info<< nl
1032  << "Disabled writing of merged mesh (-addressing-only)"
1033  << nl << nl;
1034  }
1035  else
1036  {
1037  // Write reconstructed mesh
1038  writeMesh(masterMesh[0], cellProcAddressing);
1039  if (writeCellDist)
1040  {
1041  writeDistribution
1042  (
1043  runTime,
1044  masterMesh[0],
1045  cellProcAddressing
1046  );
1047  }
1048  }
1049  }
1050  else
1051  {
1052  // Load all meshes
1053  PtrList<fvMesh> fvMeshes(nProcs);
1054  for (label proci=0; proci<nProcs; proci++)
1055  {
1056  fvMeshes.set
1057  (
1058  proci,
1059  new fvMesh
1060  (
1061  IOobject
1062  (
1063  regionName,
1064  databases[proci].timeName(),
1065  databases[proci]
1066  )
1067  )
1068  );
1069  }
1070 
1071  // Construct pointers to meshes
1072  UPtrList<polyMesh> meshes(fvMeshes.size());
1073  forAll(fvMeshes, proci)
1074  {
1075  meshes.set(proci, &fvMeshes[proci]);
1076  }
1077 
1078  // Get pairs of patches to stitch. These pairs have to
1079  // - have ordered, opposite faces (so one to one correspondence)
1080  List<DynamicList<label>> localPatch;
1081  List<DynamicList<label>> remoteProc;
1082  List<DynamicList<label>> remotePatch;
1083  const label nGlobalPatches = polyMeshAdder::procPatchPairs
1084  (
1085  meshes,
1086  localPatch,
1087  remoteProc,
1088  remotePatch
1089  );
1090 
1091  // Collect matching boundary faces on patches-to-stitch
1092  labelListList localBoundaryFace;
1093  labelListList remoteFaceProc;
1094  labelListList remoteBoundaryFace;
1096  (
1097  meshes,
1098  localPatch,
1099  remoteProc,
1100  remotePatch,
1101  localBoundaryFace,
1102  remoteFaceProc,
1103  remoteBoundaryFace
1104  );
1105 
1106  // All matched faces assumed to have vertex0 matched
1107  labelListList remoteFaceStart(meshes.size());
1108  forAll(meshes, proci)
1109  {
1110  const labelList& procFaces = localBoundaryFace[proci];
1111  remoteFaceStart[proci].setSize(procFaces.size(), 0);
1112  }
1113 
1114 
1115 
1116  labelListList patchMap(meshes.size());
1117  labelListList pointZoneMap(meshes.size());
1118  labelListList faceZoneMap(meshes.size());
1119  labelListList cellZoneMap(meshes.size());
1120  forAll(meshes, proci)
1121  {
1122  const polyMesh& mesh = meshes[proci];
1123  patchMap[proci] = identity(mesh.boundaryMesh().size());
1124 
1125  // Remove excess patches
1126  patchMap[proci].setSize(nGlobalPatches);
1127 
1128  pointZoneMap[proci] = identity(mesh.pointZones().size());
1129  faceZoneMap[proci] = identity(mesh.faceZones().size());
1130  cellZoneMap[proci] = identity(mesh.cellZones().size());
1131  }
1132 
1133 
1134  refPtr<fvMesh> masterMeshPtr;
1135  {
1136  // Do in-place addition on proc0.
1137 
1138  const labelList oldFaceOwner(fvMeshes[0].faceOwner());
1139 
1141  (
1142  0, // index of mesh to modify (== mesh_)
1143  fvMeshes,
1144  oldFaceOwner,
1145 
1146  // Coupling info
1147  localBoundaryFace,
1148  remoteFaceProc,
1149  remoteBoundaryFace,
1150 
1151  boundProcAddressing,
1152  cellProcAddressing,
1153  faceProcAddressing,
1154  pointProcAddressing
1155  );
1156 
1157  // Remove zero-faces processor patches
1158  const polyBoundaryMesh& pbm = fvMeshes[0].boundaryMesh();
1159  labelList oldToNew(pbm.size(), -1);
1160  label newi = 0;
1161  // Non processor patches first
1162  forAll(pbm, patchi)
1163  {
1164  const auto& pp = pbm[patchi];
1165  if (!isA<processorPolyPatch>(pp) || pp.size())
1166  {
1167  oldToNew[patchi] = newi++;
1168  }
1169  }
1170  const label nNonProcPatches = newi;
1171 
1172  // Move all deletable patches to the end
1173  forAll(oldToNew, patchi)
1174  {
1175  if (oldToNew[patchi] == -1)
1176  {
1177  oldToNew[patchi] = newi++;
1178  }
1179  }
1181  (
1182  fvMeshes[0],
1183  oldToNew,
1184  nNonProcPatches,
1185  false
1186  );
1187 
1188  masterMeshPtr.cref(fvMeshes[0]);
1189  }
1190 
1191 
1192  const fvMesh& masterMesh = masterMeshPtr();
1193 
1194  // Number of internal faces on the final reconstructed mesh
1195  masterInternalFaces = masterMesh.nInternalFaces();
1196 
1197  // Owner addressing on the final reconstructed mesh
1198  masterOwner = masterMesh.faceOwner();
1199 
1200 
1201  // Write reconstructed mesh
1202  // Override:
1203  // - caseName
1204  // - processorCase flag
1205  // so the resulting mesh goes to the correct location (even with
1206  // collated). The better way of solving this is to construct
1207  // (zero) mesh on the undecomposed runTime.
1208 
1209  if (writeAddrOnly)
1210  {
1211  Info<< nl
1212  << "Disabled writing of merged mesh (-addressing-only)"
1213  << nl << nl;
1214  }
1215  else
1216  {
1217  Time& masterTime = const_cast<Time&>(masterMesh.time());
1218 
1219  const word oldCaseName = masterTime.caseName();
1220  masterTime.caseName() = runTime.caseName();
1221  const bool oldProcCase(masterTime.processorCase(false));
1222 
1223  // Write reconstructed mesh
1224  writeMesh(masterMesh, cellProcAddressing);
1225  if (writeCellDist)
1226  {
1227  writeDistribution
1228  (
1229  runTime,
1230  masterMesh,
1231  cellProcAddressing
1232  );
1233  }
1234  masterTime.caseName() = oldCaseName;
1235  masterTime.processorCase(oldProcCase);
1236  }
1237  }
1238 
1239 
1240  // Write the addressing
1241 
1242  Info<< "Reconstructing addressing from processor meshes"
1243  << " to the newly reconstructed mesh" << nl << endl;
1244 
1245  forAll(databases, proci)
1246  {
1247  Info<< "Processor " << proci << nl
1248  << "Read processor mesh: "
1249  << (databases[proci].caseName()/regionDir) << endl;
1250 
1251  polyMesh procMesh
1252  (
1253  IOobject
1254  (
1255  regionName,
1256  databases[proci].timeName(),
1257  databases[proci]
1258  )
1259  );
1260 
1261  writeMaps
1262  (
1263  masterInternalFaces,
1264  masterOwner,
1265  procMesh,
1266  cellProcAddressing[proci],
1267  faceProcAddressing[proci],
1268  pointProcAddressing[proci],
1269  boundProcAddressing[proci]
1270  );
1271  }
1272  }
1273  }
1274 
1275  Info<< "End\n" << endl;
1276 
1277  return 0;
1278 }
1279 
1280 
1281 // ************************************************************************* //
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:598
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 (not the indices) of 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:666
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
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:423
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:2088
#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:2106
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:618
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:316
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:994
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:608
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:405
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:935
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:1102
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:670
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:662
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:678
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:74
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:172
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