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 | Copyright (C) 2011-2017 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 
54 using namespace Foam;
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 // Tolerance (as fraction of the bounding box). Needs to be fairly lax since
59 // usually meshes get written with limited precision (6 digits)
60 static const scalar defaultMergeTol = 1e-7;
61 
62 
63 static void renumber
64 (
65  const labelList& map,
66  labelList& elems
67 )
68 {
69  forAll(elems, i)
70  {
71  if (elems[i] >= 0)
72  {
73  elems[i] = map[elems[i]];
74  }
75  }
76 }
77 
78 
79 // Determine which faces are coupled. Uses geometric merge distance.
80 // Looks either at all boundaryFaces (fullMatch) or only at the
81 // procBoundaries for proci. Assumes that masterMesh contains already merged
82 // all the processors < proci.
83 autoPtr<faceCoupleInfo> determineCoupledFaces
84 (
85  const bool fullMatch,
86  const label masterMeshProcStart,
87  const label masterMeshProcEnd,
88  const polyMesh& masterMesh,
89  const label meshToAddProcStart,
90  const label meshToAddProcEnd,
91  const polyMesh& meshToAdd,
92  const scalar mergeDist
93 )
94 {
95  if (fullMatch || masterMesh.nCells() == 0)
96  {
98  (
99  new faceCoupleInfo
100  (
101  masterMesh,
102  meshToAdd,
103  mergeDist, // Absolute merging distance
104  true // Matching faces identical
105  )
106  );
107  }
108  else
109  {
110  // Pick up all patches on masterMesh ending in "toDDD" where DDD is
111  // the processor number proci.
112 
113  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
114 
115 
116  DynamicList<label> masterFaces
117  (
118  masterMesh.nFaces()
119  - masterMesh.nInternalFaces()
120  );
121 
122 
123  forAll(masterPatches, patchi)
124  {
125  const polyPatch& pp = masterPatches[patchi];
126 
127  if (isA<processorPolyPatch>(pp))
128  {
129  for
130  (
131  label proci=meshToAddProcStart;
132  proci<meshToAddProcEnd;
133  proci++
134  )
135  {
136  const string toProcString("to" + name(proci));
137  if (
138  pp.name().rfind(toProcString)
139  == (pp.name().size()-toProcString.size())
140  )
141  {
142  label meshFacei = pp.start();
143  forAll(pp, i)
144  {
145  masterFaces.append(meshFacei++);
146  }
147  break;
148  }
149  }
150 
151  }
152  }
153  masterFaces.shrink();
154 
155 
156  // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
157  // where DDD is the processor number proci and YYY is < proci.
158 
159  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
160 
161  DynamicList<label> addFaces
162  (
163  meshToAdd.nFaces()
164  - meshToAdd.nInternalFaces()
165  );
166 
167  forAll(addPatches, patchi)
168  {
169  const polyPatch& pp = addPatches[patchi];
170 
171  if (isA<processorPolyPatch>(pp))
172  {
173  bool isConnected = false;
174 
175  for
176  (
177  label mergedProci=masterMeshProcStart;
178  !isConnected && (mergedProci < masterMeshProcEnd);
179  mergedProci++
180  )
181  {
182  for
183  (
184  label proci = meshToAddProcStart;
185  proci < meshToAddProcEnd;
186  proci++
187  )
188  {
189  const word fromProcString
190  (
191  processorPolyPatch::newName(proci, mergedProci)
192  );
193 
194  if (pp.name() == fromProcString)
195  {
196  isConnected = true;
197  break;
198  }
199  }
200  }
201 
202  if (isConnected)
203  {
204  label meshFacei = pp.start();
205  forAll(pp, i)
206  {
207  addFaces.append(meshFacei++);
208  }
209  }
210  }
211  }
212  addFaces.shrink();
213 
215  (
216  new faceCoupleInfo
217  (
218  masterMesh,
219  masterFaces,
220  meshToAdd,
221  addFaces,
222  mergeDist, // Absolute merging distance
223  true, // Matching faces identical?
224  false, // If perfect match are faces already ordered
225  // (e.g. processor patches)
226  false // are faces each on separate patch?
227  )
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 (returnReduce(pointToMaster.size(), sumOp<label>()) == 0)
255  {
256  return autoPtr<mapPolyMesh>(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 argList& args,
308  const PtrList<Time>& databases,
309  const word& regionDir
310 )
311 {
313 
314  forAll(databases, proci)
315  {
316  fileName pointsInstance
317  (
318  databases[proci].findInstance
319  (
320  regionDir/polyMesh::meshSubDir,
321  "points"
322  )
323  );
324 
325  if (pointsInstance != databases[proci].timeName())
326  {
328  << "Your time was specified as " << databases[proci].timeName()
329  << " but there is no polyMesh/points in that time." << endl
330  << "(there is a points file in " << pointsInstance
331  << ")" << endl
332  << "Please rerun with the correct time specified"
333  << " (through the -constant, -time or -latestTime "
334  << "(at your option)."
335  << endl << exit(FatalError);
336  }
337 
338  Info<< "Reading points from "
339  << databases[proci].caseName()
340  << " for time = " << databases[proci].timeName()
341  << nl << endl;
342 
344  (
345  IOobject
346  (
347  "points",
348  databases[proci].findInstance
349  (
350  regionDir/polyMesh::meshSubDir,
351  "points"
352  ),
353  regionDir/polyMesh::meshSubDir,
354  databases[proci],
356  IOobject::NO_WRITE,
357  false
358  )
359  );
360 
361  boundBox domainBb(points, false);
362 
363  bb.min() = min(bb.min(), domainBb.min());
364  bb.max() = max(bb.max(), domainBb.max());
365  }
366 
367  return bb;
368 }
369 
370 
371 void writeCellDistance
372 (
373  Time& runTime,
374  const fvMesh& masterMesh,
375  const labelListList& cellProcAddressing
376 
377 )
378 {
379  // Write the decomposition as labelList for use with 'manual'
380  // decomposition method.
381  labelIOList cellDecomposition
382  (
383  IOobject
384  (
385  "cellDecomposition",
386  masterMesh.facesInstance(),
387  masterMesh,
390  false
391  ),
392  masterMesh.nCells()
393  );
394 
395  forAll(cellProcAddressing, proci)
396  {
397  const labelList& pCells = cellProcAddressing[proci];
398  UIndirectList<label>(cellDecomposition, pCells) = proci;
399  }
400 
401  cellDecomposition.write();
402 
403  Info<< nl << "Wrote decomposition to "
404  << cellDecomposition.objectPath()
405  << " for use in manual decomposition." << endl;
406 
407 
408  // Write as volScalarField for postprocessing. Change time to 0
409  // if was 'constant'
410  {
411  const scalar oldTime = runTime.value();
412  const label oldIndex = runTime.timeIndex();
413  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
414  {
415  runTime.setTime(0, oldIndex+1);
416  }
417 
418  volScalarField cellDist
419  (
420  IOobject
421  (
422  "cellDist",
423  runTime.timeName(),
424  masterMesh,
427  ),
428  masterMesh,
429  dimensionedScalar("cellDist", dimless, 0),
430  extrapolatedCalculatedFvPatchScalarField::typeName
431  );
432 
433  forAll(cellDecomposition, celli)
434  {
435  cellDist[celli] = cellDecomposition[celli];
436  }
437  cellDist.correctBoundaryConditions();
438 
439  cellDist.write();
440 
441  Info<< nl << "Wrote decomposition as volScalarField to "
442  << cellDist.name() << " for use in postprocessing."
443  << endl;
444 
445  // Restore time
446  runTime.setTime(oldTime, oldIndex);
447  }
448 }
449 
450 
451 int main(int argc, char *argv[])
452 {
454  (
455  "reconstruct a mesh using geometric information only"
456  );
457 
458  // Enable -constant ... if someone really wants it
459  // Enable -withZero to prevent accidentally trashing the initial fields
460  timeSelector::addOptions(true, true);
463  (
464  "mergeTol",
465  "scalar",
466  "specify the merge distance relative to the bounding box size "
467  "(default 1e-7)"
468  );
470  (
471  "fullMatch",
472  "do (slower) geometric matching on all boundary faces"
473  );
475  (
476  "cellDist",
477  "write cell distribution as a labelList - for use with 'manual' "
478  "decomposition method or as a volScalarField for post-processing."
479  );
480 
481  #include "addRegionOption.H"
482  #include "setRootCase.H"
483  #include "createTime.H"
484 
485  Info<< "This is an experimental tool which tries to merge"
486  << " individual processor" << nl
487  << "meshes back into one master mesh. Use it if the original"
488  << " master mesh has" << nl
489  << "been deleted or if the processor meshes have been modified"
490  << " (topology change)." << nl
491  << "This tool will write the resulting mesh to a new time step"
492  << " and construct" << nl
493  << "xxxxProcAddressing files in the processor meshes so"
494  << " reconstructPar can be" << nl
495  << "used to regenerate the fields on the master mesh." << nl
496  << nl
497  << "Not well tested & use at your own risk!" << nl
498  << endl;
499 
500 
501  word regionName = polyMesh::defaultRegion;
502  word regionDir = word::null;
503 
504  if
505  (
506  args.optionReadIfPresent("region", regionName)
507  && regionName != polyMesh::defaultRegion
508  )
509  {
510  regionDir = regionName;
511  Info<< "Operating on region " << regionName << nl << endl;
512  }
513 
514  scalar mergeTol = defaultMergeTol;
515  args.optionReadIfPresent("mergeTol", mergeTol);
516 
517  scalar writeTol = Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
518 
519  Info<< "Merge tolerance : " << mergeTol << nl
520  << "Write tolerance : " << writeTol << endl;
521 
522  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
523  {
525  << "Your current settings specify ASCII writing with "
526  << IOstream::defaultPrecision() << " digits precision." << endl
527  << "Your merging tolerance (" << mergeTol << ") is finer than this."
528  << endl
529  << "Please change your writeFormat to binary"
530  << " or increase the writePrecision" << endl
531  << "or adjust the merge tolerance (-mergeTol)."
532  << exit(FatalError);
533  }
534 
535 
536  const bool fullMatch = args.optionFound("fullMatch");
537 
538  if (fullMatch)
539  {
540  Info<< "Doing geometric matching on all boundary faces." << nl << endl;
541  }
542  else
543  {
544  Info<< "Doing geometric matching on correct procBoundaries only."
545  << nl << "This assumes a correct decomposition." << endl;
546  }
547 
548  bool writeCellDist = args.optionFound("cellDist");
549 
550 
551  label nProcs = fileHandler().nProcs(args.path());
552 
553  Info<< "Found " << nProcs << " processor directories" << nl << endl;
554 
555 
556  // Read all time databases
557  PtrList<Time> databases(nProcs);
558 
559  forAll(databases, proci)
560  {
561  Info<< "Reading database "
562  << args.caseName()/fileName(word("processor") + name(proci))
563  << endl;
564 
565  databases.set
566  (
567  proci,
568  new Time
569  (
571  args.rootPath(),
572  args.caseName()/fileName(word("processor") + name(proci))
573  )
574  );
575  }
576 
577  // Use the times list from the master processor
578  // and select a subset based on the command-line options
580  (
581  databases[0].times(),
582  args
583  );
584 
585  // Loop over all times
586  forAll(timeDirs, timeI)
587  {
588  // Set time for global database
589  runTime.setTime(timeDirs[timeI], timeI);
590 
591  Info<< "Time = " << runTime.timeName() << nl << endl;
592 
593  // Set time for all databases
594  forAll(databases, proci)
595  {
596  databases[proci].setTime(timeDirs[timeI], timeI);
597  }
598 
599  IOobject facesIO
600  (
601  "faces",
602  databases[0].timeName(),
603  regionDir/polyMesh::meshSubDir,
604  databases[0],
605  IOobject::NO_READ,
606  IOobject::NO_WRITE
607  );
608 
609 
610  // Problem: faceCompactIOList recognises both 'faceList' and
611  // 'faceCompactList' so we should be lenient when doing
612  // typeHeaderOk
613  if (!facesIO.typeHeaderOk<faceCompactIOList>(false))
614  {
615  Info<< "No mesh." << nl << endl;
616  continue;
617  }
618 
619 
620  // Read point on individual processors to determine merge tolerance
621  // (otherwise single cell domains might give problems)
622 
623  const boundBox bb = procBounds(args, databases, regionDir);
624  const scalar mergeDist = mergeTol*bb.mag();
625 
626  Info<< "Overall mesh bounding box : " << bb << nl
627  << "Relative tolerance : " << mergeTol << nl
628  << "Absolute matching distance : " << mergeDist << nl
629  << endl;
630 
631 
632  // Addressing from processor to reconstructed case
633  labelListList cellProcAddressing(nProcs);
635  labelListList pointProcAddressing(nProcs);
636  labelListList boundaryProcAddressing(nProcs);
637 
638  // Internal faces on the final reconstructed mesh
639  label masterInternalFaces;
640 
641  // Owner addressing on the final reconstructed mesh
642  labelList masterOwner;
643 
644  {
645  // Construct empty mesh.
646  // fvMesh** masterMesh = new fvMesh*[nProcs];
647  PtrList<fvMesh> masterMesh(nProcs);
648 
649  for (label proci=0; proci<nProcs; proci++)
650  {
651  masterMesh.set
652  (
653  proci,
654  new fvMesh
655  (
656  IOobject
657  (
658  regionName,
659  runTime.timeName(),
660  runTime,
662  ),
663  xferCopy(pointField()),
664  xferCopy(faceList()),
665  xferCopy(cellList())
666  )
667  );
668 
669  fvMesh meshToAdd
670  (
671  IOobject
672  (
673  regionName,
674  databases[proci].timeName(),
675  databases[proci]
676  )
677  );
678 
679  // Initialize its addressing
680  cellProcAddressing[proci] = identity(meshToAdd.nCells());
681  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
682  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
683  boundaryProcAddressing[proci] =
684  identity(meshToAdd.boundaryMesh().size());
685 
686  // Find geometrically shared points/faces.
687  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
688  (
689  fullMatch,
690  proci,
691  proci,
692  masterMesh[proci],
693  proci,
694  proci,
695  meshToAdd,
696  mergeDist
697  );
698 
699  // Add elements to mesh
701  (
702  masterMesh[proci],
703  meshToAdd,
704  couples
705  );
706 
707  // Added processor
708  renumber(map().addedCellMap(), cellProcAddressing[proci]);
709  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
710  renumber(map().addedPointMap(), pointProcAddressing[proci]);
711  renumber(map().addedPatchMap(), boundaryProcAddressing[proci]);
712  }
713  for (label step=2; step<nProcs*2; step*=2)
714  {
715  for (label proci=0; proci<nProcs; proci+=step)
716  {
717  label next = proci + step/2;
718  if(next >= nProcs)
719  {
720  continue;
721  }
722 
723  Info<< "Merging mesh " << proci << " with " << next << endl;
724 
725  // Find geometrically shared points/faces.
726  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
727  (
728  fullMatch,
729  proci,
730  next,
731  masterMesh[proci],
732  next,
733  proci+step,
734  masterMesh[next],
735  mergeDist
736  );
737 
738  // Add elements to mesh
740  (
741  masterMesh[proci],
742  masterMesh[next],
743  couples
744  );
745 
746  // Processors that were already in masterMesh
747  for (label mergedI=proci; mergedI<next; mergedI++)
748  {
749  renumber
750  (
751  map().oldCellMap(),
752  cellProcAddressing[mergedI]
753  );
754 
755  renumber
756  (
757  map().oldFaceMap(),
758  faceProcAddressing[mergedI]
759  );
760 
761  renumber
762  (
763  map().oldPointMap(),
764  pointProcAddressing[mergedI]
765  );
766 
767  // Note: boundary is special since can contain -1.
768  renumber
769  (
770  map().oldPatchMap(),
771  boundaryProcAddressing[mergedI]
772  );
773  }
774 
775  // Added processor
776  for
777  (
778  label addedI=next;
779  addedI<min(proci+step, nProcs);
780  addedI++
781  )
782  {
783  renumber
784  (
785  map().addedCellMap(),
786  cellProcAddressing[addedI]
787  );
788 
789  renumber
790  (
791  map().addedFaceMap(),
792  faceProcAddressing[addedI]
793  );
794 
795  renumber
796  (
797  map().addedPointMap(),
798  pointProcAddressing[addedI]
799  );
800 
801  renumber
802  (
803  map().addedPatchMap(),
804  boundaryProcAddressing[addedI]
805  );
806  }
807 
808  masterMesh.set(next, nullptr);
809  }
810  }
811 
812  for (label proci=0; proci<nProcs; proci++)
813  {
814  Info<< "Reading mesh to add from "
815  << databases[proci].caseName()
816  << " for time = " << databases[proci].timeName()
817  << nl << nl << endl;
818  }
819 
820  // See if any points on the mastermesh have become connected
821  // because of connections through processor meshes.
822  mergeSharedPoints(mergeDist, masterMesh[0], pointProcAddressing);
823 
824  // Save some properties on the reconstructed mesh
825  masterInternalFaces = masterMesh[0].nInternalFaces();
826  masterOwner = masterMesh[0].faceOwner();
827 
828 
829  Info<< "\nWriting merged mesh to "
830  << runTime.path()/runTime.timeName()
831  << nl << endl;
832 
833  if (!masterMesh[0].write())
834  {
836  << "Failed writing polyMesh."
837  << exit(FatalError);
838  }
839 
840  if (writeCellDist)
841  {
842  writeCellDistance
843  (
844  runTime,
845  masterMesh[0],
846  cellProcAddressing
847  );
848  }
849  }
850 
851 
852  // Write the addressing
853 
854  Info<< "Reconstructing the addressing from the processor meshes"
855  << " to the newly reconstructed mesh" << nl << endl;
856 
857  forAll(databases, proci)
858  {
859  Info<< "Reading processor " << proci << " mesh from "
860  << databases[proci].caseName() << endl;
861 
862  polyMesh procMesh
863  (
864  IOobject
865  (
866  regionName,
867  databases[proci].timeName(),
868  databases[proci]
869  )
870  );
871 
872 
873  // From processor point to reconstructed mesh point
874 
875  Info<< "Writing pointProcAddressing to "
876  << databases[proci].caseName()
877  /procMesh.facesInstance()
879  << endl;
880 
882  (
883  IOobject
884  (
885  "pointProcAddressing",
886  procMesh.facesInstance(),
888  procMesh,
891  false // Do not register
892  ),
893  pointProcAddressing[proci]
894  ).write();
895 
896 
897  // From processor face to reconstructed mesh face
898 
899  Info<< "Writing faceProcAddressing to "
900  << databases[proci].caseName()
901  /procMesh.facesInstance()
903  << endl;
904 
905  labelIOList faceProcAddr
906  (
907  IOobject
908  (
909  "faceProcAddressing",
910  procMesh.facesInstance(),
912  procMesh,
915  false // Do not register
916  ),
917  faceProcAddressing[proci]
918  );
919 
920  // Now add turning index to faceProcAddressing.
921  // See reconstructPar for meaning of turning index.
922  forAll(faceProcAddr, procFacei)
923  {
924  const label masterFacei = faceProcAddr[procFacei];
925 
926  if
927  (
928  !procMesh.isInternalFace(procFacei)
929  && masterFacei < masterInternalFaces
930  )
931  {
932  // proc face is now external but used to be internal face.
933  // Check if we have owner or neighbour.
934 
935  label procOwn = procMesh.faceOwner()[procFacei];
936  label masterOwn = masterOwner[masterFacei];
937 
938  if (cellProcAddressing[proci][procOwn] == masterOwn)
939  {
940  // No turning. Offset by 1.
941  faceProcAddr[procFacei]++;
942  }
943  else
944  {
945  // Turned face.
946  faceProcAddr[procFacei] =
947  -1 - faceProcAddr[procFacei];
948  }
949  }
950  else
951  {
952  // No turning. Offset by 1.
953  faceProcAddr[procFacei]++;
954  }
955  }
956 
957  faceProcAddr.write();
958 
959 
960  // From processor cell to reconstructed mesh cell
961 
962  Info<< "Writing cellProcAddressing to "
963  << databases[proci].caseName()
964  /procMesh.facesInstance()
966  << endl;
967 
969  (
970  IOobject
971  (
972  "cellProcAddressing",
973  procMesh.facesInstance(),
975  procMesh,
978  false // Do not register
979  ),
980  cellProcAddressing[proci]
981  ).write();
982 
983 
984 
985  // From processor patch to reconstructed mesh patch
986 
987  Info<< "Writing boundaryProcAddressing to "
988  << databases[proci].caseName()
989  /procMesh.facesInstance()
991  << endl;
992 
994  (
995  IOobject
996  (
997  "boundaryProcAddressing",
998  procMesh.facesInstance(),
1000  procMesh,
1003  false // Do not register
1004  ),
1005  boundaryProcAddressing[proci]
1006  ).write();
1007 
1008  Info<< endl;
1009  }
1010  }
1011 
1012 
1013  Info<< "End.\n" << endl;
1014 
1015  return 0;
1016 }
1017 
1018 
1019 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
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:428
fileName path() const
Return path.
Definition: Time.H:269
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:69
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
static const Form max
Definition: VectorSpace.H:112
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:801
const double e
Elementary charge.
Definition: doubleFloat.H:78
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
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:253
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:148
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
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
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:113
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:95
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
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:856
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:201
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:262
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:291
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:63
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:300
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
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:85
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:126
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual bool write(const bool valid=true) const
Write using setting from DB.
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
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.