reconstructParMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  reconstructParMesh
26 
27 Description
28  Reconstructs a mesh using geometric information only.
29 
30  Writes point/face/cell procAddressing so afterwards reconstructPar can be
31  used to reconstruct fields.
32 
33  Note:
34  - uses geometric matching tolerance (set with -mergeTol (at your option)
35 
36  If the parallel case does not have correct procBoundaries use the
37  -fullMatch option which will check all boundary faces (bit slower).
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "argList.H"
42 #include "timeSelector.H"
43 
44 #include "IOobjectList.H"
45 #include "labelIOList.H"
46 #include "processorPolyPatch.H"
47 #include "mapAddedPolyMesh.H"
48 #include "polyMeshAdder.H"
49 #include "faceCoupleInfo.H"
50 #include "fvMeshAdder.H"
51 #include "polyTopoChange.H"
53 #include "regionProperties.H"
54 
55 using namespace Foam;
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
60 // usually meshes get written with limited precision (6 digits)
61 static const scalar defaultMergeTol = 1e-7;
62 
63 
64 static void renumber
65 (
66  const labelList& map,
67  labelList& elems
68 )
69 {
70  forAll(elems, i)
71  {
72  if (elems[i] >= 0)
73  {
74  elems[i] = map[elems[i]];
75  }
76  }
77 }
78 
79 
80 // Determine which faces are coupled. Uses geometric merge distance.
81 // Looks either at all boundaryFaces (fullMatch) or only at the
82 // procBoundaries for proci. Assumes that masterMesh contains already merged
83 // all the processors < proci.
84 autoPtr<faceCoupleInfo> determineCoupledFaces
85 (
86  const bool fullMatch,
87  const label masterMeshProcStart,
88  const label masterMeshProcEnd,
89  const polyMesh& masterMesh,
90  const label meshToAddProcStart,
91  const label meshToAddProcEnd,
92  const polyMesh& meshToAdd,
93  const scalar mergeDist
94 )
95 {
96  if (fullMatch || masterMesh.nCells() == 0)
97  {
99  (
100  new faceCoupleInfo
101  (
102  masterMesh,
103  meshToAdd,
104  mergeDist, // Absolute merging distance
105  true // Matching faces identical
106  )
107  );
108  }
109  else
110  {
111  // Pick up all patches on masterMesh ending in "toDDD" where DDD is
112  // the processor number proci.
113 
114  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
115 
116 
117  DynamicList<label> masterFaces
118  (
119  masterMesh.nFaces()
120  - masterMesh.nInternalFaces()
121  );
122 
123 
124  forAll(masterPatches, patchi)
125  {
126  const polyPatch& pp = masterPatches[patchi];
127 
128  if (isA<processorPolyPatch>(pp))
129  {
130  for
131  (
132  label proci=meshToAddProcStart;
133  proci<meshToAddProcEnd;
134  proci++
135  )
136  {
137  const string toProcString("to" + name(proci));
138  if (
139  pp.name().rfind(toProcString)
140  == (pp.name().size()-toProcString.size())
141  )
142  {
143  label meshFacei = pp.start();
144  forAll(pp, i)
145  {
146  masterFaces.append(meshFacei++);
147  }
148  break;
149  }
150  }
151 
152  }
153  }
154  masterFaces.shrink();
155 
156 
157  // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
158  // where DDD is the processor number proci and YYY is < proci.
159 
160  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
161 
162  DynamicList<label> addFaces
163  (
164  meshToAdd.nFaces()
165  - meshToAdd.nInternalFaces()
166  );
167 
168  forAll(addPatches, patchi)
169  {
170  const polyPatch& pp = addPatches[patchi];
171 
172  if (isA<processorPolyPatch>(pp))
173  {
174  bool isConnected = false;
175 
176  for
177  (
178  label mergedProci=masterMeshProcStart;
179  !isConnected && (mergedProci < masterMeshProcEnd);
180  mergedProci++
181  )
182  {
183  for
184  (
185  label proci = meshToAddProcStart;
186  proci < meshToAddProcEnd;
187  proci++
188  )
189  {
190  const word fromProcString
191  (
192  processorPolyPatch::newName(proci, mergedProci)
193  );
194 
195  if (pp.name() == fromProcString)
196  {
197  isConnected = true;
198  break;
199  }
200  }
201  }
202 
203  if (isConnected)
204  {
205  label meshFacei = pp.start();
206  forAll(pp, i)
207  {
208  addFaces.append(meshFacei++);
209  }
210  }
211  }
212  }
213  addFaces.shrink();
214 
216  (
217  new faceCoupleInfo
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 
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 autoPtr<mapPolyMesh>(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 argList& args,
309  const PtrList<Time>& databases,
310  const word& regionDir
311 )
312 {
314 
315  forAll(databases, proci)
316  {
317  fileName pointsInstance
318  (
319  databases[proci].findInstance
320  (
321  regionDir/polyMesh::meshSubDir,
322  "points"
323  )
324  );
325 
326  if (pointsInstance != databases[proci].timeName())
327  {
329  << "Your time was specified as " << databases[proci].timeName()
330  << " but there is no polyMesh/points in that time." << endl
331  << "(there is a points file in " << pointsInstance
332  << ")" << endl
333  << "Please rerun with the correct time specified"
334  << " (through the -constant, -time or -latestTime "
335  << "(at your option)."
336  << endl << exit(FatalError);
337  }
338 
339  Info<< "Reading points from "
340  << databases[proci].caseName()
341  << " for time = " << databases[proci].timeName()
342  << nl << endl;
343 
345  (
346  IOobject
347  (
348  "points",
349  databases[proci].findInstance
350  (
351  regionDir/polyMesh::meshSubDir,
352  "points"
353  ),
354  regionDir/polyMesh::meshSubDir,
355  databases[proci],
357  IOobject::NO_WRITE,
358  false
359  )
360  );
361 
362  boundBox domainBb(points, false);
363 
364  bb.min() = min(bb.min(), domainBb.min());
365  bb.max() = max(bb.max(), domainBb.max());
366  }
367 
368  return bb;
369 }
370 
371 
372 void writeCellDistance
373 (
374  Time& runTime,
375  const fvMesh& masterMesh,
376  const labelListList& cellProcAddressing
377 
378 )
379 {
380  // Write the decomposition as labelList for use with 'manual'
381  // decomposition method.
382  labelIOList cellDecomposition
383  (
384  IOobject
385  (
386  "cellDecomposition",
387  masterMesh.facesInstance(),
388  masterMesh,
391  false
392  ),
393  masterMesh.nCells()
394  );
395 
396  forAll(cellProcAddressing, proci)
397  {
398  const labelList& pCells = cellProcAddressing[proci];
399  UIndirectList<label>(cellDecomposition, pCells) = proci;
400  }
401 
402  cellDecomposition.write();
403 
404  Info<< nl << "Wrote decomposition to "
405  << cellDecomposition.objectPath()
406  << " for use in manual decomposition." << endl;
407 
408 
409  // Write as volScalarField for postprocessing. Change time to 0
410  // if was 'constant'
411  {
412  const scalar oldTime = runTime.value();
413  const label oldIndex = runTime.timeIndex();
414  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
415  {
416  runTime.setTime(0, oldIndex+1);
417  }
418 
419  volScalarField cellDist
420  (
421  IOobject
422  (
423  "cellDist",
424  runTime.timeName(),
425  masterMesh,
428  ),
429  masterMesh,
431  extrapolatedCalculatedFvPatchScalarField::typeName
432  );
433 
434  forAll(cellDecomposition, celli)
435  {
436  cellDist[celli] = cellDecomposition[celli];
437  }
438  cellDist.correctBoundaryConditions();
439 
440  cellDist.write();
441 
442  Info<< nl << "Wrote decomposition as volScalarField to "
443  << cellDist.name() << " for use in postprocessing."
444  << endl;
445 
446  // Restore time
447  runTime.setTime(oldTime, oldIndex);
448  }
449 }
450 
451 
452 int main(int argc, char *argv[])
453 {
455  (
456  "reconstruct a mesh using geometric information only"
457  );
458 
459  // Enable -constant ... if someone really wants it
460  // Enable -withZero to prevent accidentally trashing the initial fields
461  timeSelector::addOptions(true, true);
464  (
465  "mergeTol",
466  "scalar",
467  "specify the merge distance relative to the bounding box size "
468  "(default 1e-7)"
469  );
471  (
472  "fullMatch",
473  "do (slower) geometric matching on all boundary faces"
474  );
476  (
477  "cellDist",
478  "write cell distribution as a labelList - for use with 'manual' "
479  "decomposition method or as a volScalarField for post-processing."
480  );
481 
482  #include "addRegionOption.H"
483  #include "addAllRegionsOption.H"
484  #include "setRootCase.H"
485  #include "createTime.H"
486 
487  Info<< "This is an experimental tool which tries to merge"
488  << " individual processor" << nl
489  << "meshes back into one master mesh. Use it if the original"
490  << " master mesh has" << nl
491  << "been deleted or if the processor meshes have been modified"
492  << " (topology change)." << nl
493  << "This tool will write the resulting mesh to a new time step"
494  << " and construct" << nl
495  << "xxxxProcAddressing files in the processor meshes so"
496  << " reconstructPar can be" << nl
497  << "used to regenerate the fields on the master mesh." << nl
498  << nl
499  << "Not well tested & use at your own risk!" << nl
500  << endl;
501 
502 
503  const wordList regionNames(selectRegionNames(args, runTime));
504  if (regionNames.size() > 1)
505  {
506  Info<< "Operating on regions " << regionNames[0];
507  for (label regioni = 1; regioni < regionNames.size() - 1; ++ regioni)
508  {
509  Info<< ", " << regionNames[regioni];
510  }
511  Info<< " and " << regionNames.last() << nl << endl;
512  }
513  else if (regionNames[0] != polyMesh::defaultRegion)
514  {
515  Info<< "Operating on region " << regionNames[0] << nl << endl;
516  }
517 
518 
519  scalar mergeTol = defaultMergeTol;
520  args.optionReadIfPresent("mergeTol", mergeTol);
521 
522  scalar writeTol = Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
523 
524  Info<< "Merge tolerance : " << mergeTol << nl
525  << "Write tolerance : " << writeTol << endl;
526 
527  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
528  {
530  << "Your current settings specify ASCII writing with "
531  << IOstream::defaultPrecision() << " digits precision." << endl
532  << "Your merging tolerance (" << mergeTol << ") is finer than this."
533  << endl
534  << "Please change your writeFormat to binary"
535  << " or increase the writePrecision" << endl
536  << "or adjust the merge tolerance (-mergeTol)."
537  << exit(FatalError);
538  }
539 
540 
541  const bool fullMatch = args.optionFound("fullMatch");
542 
543  if (fullMatch)
544  {
545  Info<< "Doing geometric matching on all boundary faces." << nl << endl;
546  }
547  else
548  {
549  Info<< "Doing geometric matching on correct procBoundaries only."
550  << nl << "This assumes a correct decomposition." << endl;
551  }
552 
553  bool writeCellDist = args.optionFound("cellDist");
554 
555 
556  label nProcs = fileHandler().nProcs(args.path());
557 
558  Info<< "Found " << nProcs << " processor directories" << nl << endl;
559 
560 
561  // Read all time databases
562  PtrList<Time> databases(nProcs);
563 
564  forAll(databases, proci)
565  {
566  Info<< "Reading database "
567  << args.caseName()/fileName(word("processor") + name(proci))
568  << endl;
569 
570  databases.set
571  (
572  proci,
573  new Time
574  (
576  args.rootPath(),
577  args.caseName()/fileName(word("processor") + name(proci))
578  )
579  );
580  }
581 
582  // Use the times list from the master processor
583  // and select a subset based on the command-line options
585  (
586  databases[0].times(),
587  args
588  );
589 
590  // Loop over all times
591  forAll(timeDirs, timeI)
592  {
593  // Set time for global database
594  runTime.setTime(timeDirs[timeI], timeI);
595 
596  Info<< "Time = " << runTime.timeName() << nl << endl;
597 
598  // Set time for all databases
599  forAll(databases, proci)
600  {
601  databases[proci].setTime(timeDirs[timeI], timeI);
602  }
603 
604  forAll(regionNames, regioni)
605  {
606  const word& regionName = regionNames[regioni];
607  const word regionDir =
608  regionName == polyMesh::defaultRegion
609  ? word::null
610  : regionName;
611 
612  IOobject facesIO
613  (
614  "faces",
615  databases[0].timeName(),
616  regionDir/polyMesh::meshSubDir,
617  databases[0],
618  IOobject::NO_READ,
619  IOobject::NO_WRITE
620  );
621 
622 
623  // Problem: faceCompactIOList recognises both 'faceList' and
624  // 'faceCompactList' so we should be lenient when doing
625  // typeHeaderOk
626  if (!facesIO.typeHeaderOk<faceCompactIOList>(false))
627  {
628  Info<< "No mesh." << nl << endl;
629  continue;
630  }
631 
632 
633  // Read point on individual processors to determine merge tolerance
634  // (otherwise single cell domains might give problems)
635 
636  const boundBox bb = procBounds(args, databases, regionDir);
637  const scalar mergeDist = mergeTol*bb.mag();
638 
639  Info<< "Overall mesh bounding box : " << bb << nl
640  << "Relative tolerance : " << mergeTol << nl
641  << "Absolute matching distance : " << mergeDist << nl
642  << endl;
643 
644 
645  // Addressing from processor to reconstructed case
646  labelListList cellProcAddressing(nProcs);
648  labelListList pointProcAddressing(nProcs);
649  labelListList boundaryProcAddressing(nProcs);
650 
651  // Internal faces on the final reconstructed mesh
652  label masterInternalFaces;
653 
654  // Owner addressing on the final reconstructed mesh
655  labelList masterOwner;
656 
657  {
658  // Construct empty mesh.
659  // fvMesh** masterMesh = new fvMesh*[nProcs];
660  PtrList<fvMesh> masterMesh(nProcs);
661 
662  for (label proci=0; proci<nProcs; proci++)
663  {
664  masterMesh.set
665  (
666  proci,
667  new fvMesh
668  (
669  IOobject
670  (
671  regionName,
672  runTime.timeName(),
673  runTime,
675  ),
676  pointField(),
677  faceList(),
678  cellList()
679  )
680  );
681 
682  fvMesh meshToAdd
683  (
684  IOobject
685  (
686  regionName,
687  databases[proci].timeName(),
688  databases[proci]
689  )
690  );
691 
692  // Initialize its addressing
693  cellProcAddressing[proci] = identity(meshToAdd.nCells());
694  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
695  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
696  boundaryProcAddressing[proci] =
697  identity(meshToAdd.boundaryMesh().size());
698 
699  // Find geometrically shared points/faces.
700  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
701  (
702  fullMatch,
703  proci,
704  proci,
705  masterMesh[proci],
706  proci,
707  proci,
708  meshToAdd,
709  mergeDist
710  );
711 
712  // Add elements to mesh
714  (
715  masterMesh[proci],
716  meshToAdd,
717  couples
718  );
719 
720  // Added processor
721  renumber(map().addedCellMap(), cellProcAddressing[proci]);
722  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
723  renumber(map().addedPointMap(), pointProcAddressing[proci]);
724  renumber
725  (
726  map().addedPatchMap(),
727  boundaryProcAddressing[proci]
728  );
729  }
730  for (label step=2; step<nProcs*2; step*=2)
731  {
732  for (label proci=0; proci<nProcs; proci+=step)
733  {
734  label next = proci + step/2;
735  if(next >= nProcs)
736  {
737  continue;
738  }
739 
740  Info<< "Merging mesh " << proci << " with " << next
741  << endl;
742 
743  // Find geometrically shared points/faces.
744  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
745  (
746  fullMatch,
747  proci,
748  next,
749  masterMesh[proci],
750  next,
751  proci+step,
752  masterMesh[next],
753  mergeDist
754  );
755 
756  // Add elements to mesh
758  (
759  masterMesh[proci],
760  masterMesh[next],
761  couples
762  );
763 
764  // Processors that were already in masterMesh
765  for (label mergedI=proci; mergedI<next; mergedI++)
766  {
767  renumber
768  (
769  map().oldCellMap(),
770  cellProcAddressing[mergedI]
771  );
772 
773  renumber
774  (
775  map().oldFaceMap(),
776  faceProcAddressing[mergedI]
777  );
778 
779  renumber
780  (
781  map().oldPointMap(),
782  pointProcAddressing[mergedI]
783  );
784 
785  // Note: boundary is special since can contain -1.
786  renumber
787  (
788  map().oldPatchMap(),
789  boundaryProcAddressing[mergedI]
790  );
791  }
792 
793  // Added processor
794  for
795  (
796  label addedI=next;
797  addedI<min(proci+step, nProcs);
798  addedI++
799  )
800  {
801  renumber
802  (
803  map().addedCellMap(),
804  cellProcAddressing[addedI]
805  );
806 
807  renumber
808  (
809  map().addedFaceMap(),
810  faceProcAddressing[addedI]
811  );
812 
813  renumber
814  (
815  map().addedPointMap(),
816  pointProcAddressing[addedI]
817  );
818 
819  renumber
820  (
821  map().addedPatchMap(),
822  boundaryProcAddressing[addedI]
823  );
824  }
825 
826  masterMesh.set(next, nullptr);
827  }
828  }
829 
830  for (label proci=0; proci<nProcs; proci++)
831  {
832  Info<< "Reading mesh to add from "
833  << databases[proci].caseName()
834  << " for time = " << databases[proci].timeName()
835  << nl << nl << endl;
836  }
837 
838  // See if any points on the mastermesh have become connected
839  // because of connections through processor meshes.
840  mergeSharedPoints
841  (
842  mergeDist,
843  masterMesh[0],
844  pointProcAddressing
845  );
846 
847  // Save some properties on the reconstructed mesh
848  masterInternalFaces = masterMesh[0].nInternalFaces();
849  masterOwner = masterMesh[0].faceOwner();
850 
851 
852  Info<< "\nWriting merged mesh to "
853  << runTime.path()/runTime.timeName()
854  << nl << endl;
855 
856  if (!masterMesh[0].write())
857  {
859  << "Failed writing polyMesh."
860  << exit(FatalError);
861  }
862 
863  if (writeCellDist)
864  {
865  writeCellDistance
866  (
867  runTime,
868  masterMesh[0],
869  cellProcAddressing
870  );
871  }
872  }
873 
874 
875  // Write the addressing
876 
877  Info<< "Reconstructing the addressing from the processor meshes"
878  << " to the newly reconstructed mesh" << nl << endl;
879 
880  forAll(databases, proci)
881  {
882  Info<< "Reading processor " << proci << " mesh from "
883  << databases[proci].caseName() << endl;
884 
885  polyMesh procMesh
886  (
887  IOobject
888  (
889  regionName,
890  databases[proci].timeName(),
891  databases[proci]
892  )
893  );
894 
895 
896  // From processor point to reconstructed mesh point
897 
898  Info<< "Writing pointProcAddressing to "
899  << databases[proci].caseName()
900  /procMesh.facesInstance()
902  << endl;
903 
905  (
906  IOobject
907  (
908  "pointProcAddressing",
909  procMesh.facesInstance(),
911  procMesh,
914  false // Do not register
915  ),
916  pointProcAddressing[proci]
917  ).write();
918 
919 
920  // From processor face to reconstructed mesh face
921 
922  Info<< "Writing faceProcAddressing to "
923  << databases[proci].caseName()
924  /procMesh.facesInstance()
926  << endl;
927 
928  labelIOList faceProcAddr
929  (
930  IOobject
931  (
932  "faceProcAddressing",
933  procMesh.facesInstance(),
935  procMesh,
938  false // Do not register
939  ),
940  faceProcAddressing[proci]
941  );
942 
943  // Now add turning index to faceProcAddressing.
944  // See reconstructPar for meaning of turning index.
945  forAll(faceProcAddr, procFacei)
946  {
947  const label masterFacei = faceProcAddr[procFacei];
948 
949  if
950  (
951  !procMesh.isInternalFace(procFacei)
952  && masterFacei < masterInternalFaces
953  )
954  {
955  // proc face is now external but used to be internal
956  // face. Check if we have owner or neighbour.
957 
958  label procOwn = procMesh.faceOwner()[procFacei];
959  label masterOwn = masterOwner[masterFacei];
960 
961  if (cellProcAddressing[proci][procOwn] == masterOwn)
962  {
963  // No turning. Offset by 1.
964  faceProcAddr[procFacei]++;
965  }
966  else
967  {
968  // Turned face.
969  faceProcAddr[procFacei] =
970  -1 - faceProcAddr[procFacei];
971  }
972  }
973  else
974  {
975  // No turning. Offset by 1.
976  faceProcAddr[procFacei]++;
977  }
978  }
979 
980  faceProcAddr.write();
981 
982 
983  // From processor cell to reconstructed mesh cell
984 
985  Info<< "Writing cellProcAddressing to "
986  << databases[proci].caseName()
987  /procMesh.facesInstance()
989  << endl;
990 
992  (
993  IOobject
994  (
995  "cellProcAddressing",
996  procMesh.facesInstance(),
998  procMesh,
1001  false // Do not register
1002  ),
1003  cellProcAddressing[proci]
1004  ).write();
1005 
1006 
1007 
1008  // From processor patch to reconstructed mesh patch
1009 
1010  Info<< "Writing boundaryProcAddressing to "
1011  << databases[proci].caseName()
1012  /procMesh.facesInstance()
1014  << endl;
1015 
1016  labelIOList
1017  (
1018  IOobject
1019  (
1020  "boundaryProcAddressing",
1021  procMesh.facesInstance(),
1023  procMesh,
1026  false // Do not register
1027  ),
1028  boundaryProcAddressing[proci]
1029  ).write();
1030 
1031  Info<< endl;
1032  }
1033  }
1034  }
1035 
1036 
1037  Info<< "End.\n" << endl;
1038 
1039  return 0;
1040 }
1041 
1042 
1043 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
fileName path() const
Return path.
Definition: Time.H:266
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
A class for handling file names.
Definition: fileName.H:79
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
static const Form max
Definition: VectorSpace.H:115
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:808
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nInternalFaces() const
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
PtrList< labelIOList > & faceProcAddressing
label nFaces() const
Foam::word regionName
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static void noParallel()
Remove the parallel options.
Definition: argList.C:174
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static const Form min
Definition: VectorSpace.H:116
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const pointField & points
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
A class for handling words, derived from string.
Definition: word.H:59
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:127
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const Type & value() const
Return const reference to value.
static const word null
An empty word.
Definition: word.H:77
word timeName
Definition: getTimeIndex.H:3
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:873
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:198
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
Definition: fvMeshAdder.C:71
Foam::polyBoundaryMesh.
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
static const char nl
Definition: Ostream.H:265
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:288
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName path() const
Return the path to the caseName.
Definition: argListI.H:60
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
bool set(const Key &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- vGreat.
Definition: boundBox.H:82
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
label patchi
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
label newPointi
Definition: readKivaGrid.H:501
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
messageStream Info
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
virtual bool write(const bool write=true) const
Write using setting from DB.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
List< cell > cellList
list of cells
Definition: cellList.H:42
wordList selectRegionNames(const argList &args, const Time &runTime)
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.