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