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-2015 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.
84 (
85  const bool fullMatch,
86  const label procI,
87  const polyMesh& masterMesh,
88  const polyMesh& meshToAdd,
89  const scalar mergeDist
90 )
91 {
92  if (fullMatch || masterMesh.nCells() == 0)
93  {
95  (
96  new faceCoupleInfo
97  (
98  masterMesh,
99  meshToAdd,
100  mergeDist, // Absolute merging distance
101  true // Matching faces identical
102  )
103  );
104  }
105  else
106  {
107  // Pick up all patches on masterMesh ending in "toDDD" where DDD is
108  // the processor number procI.
109 
110  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
111 
112  const string toProcString("to" + name(procI));
113 
114  DynamicList<label> masterFaces
115  (
116  masterMesh.nFaces()
117  - masterMesh.nInternalFaces()
118  );
119 
120  forAll(masterPatches, patchI)
121  {
122  const polyPatch& pp = masterPatches[patchI];
123 
124  if
125  (
126  isA<processorPolyPatch>(pp)
127  && (
128  pp.name().rfind(toProcString)
129  == (pp.name().size()-toProcString.size())
130  )
131  )
132  {
133  label meshFaceI = pp.start();
134  forAll(pp, i)
135  {
136  masterFaces.append(meshFaceI++);
137  }
138  }
139  }
140  masterFaces.shrink();
141 
142 
143  // Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
144  // where DDD is the processor number procI and YYY is < procI.
145 
146  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
147 
148  DynamicList<label> addFaces
149  (
150  meshToAdd.nFaces()
151  - meshToAdd.nInternalFaces()
152  );
153 
154  forAll(addPatches, patchI)
155  {
156  const polyPatch& pp = addPatches[patchI];
157 
158  if (isA<processorPolyPatch>(pp))
159  {
160  bool isConnected = false;
161 
162  for (label mergedProcI = 0; mergedProcI < procI; mergedProcI++)
163  {
164  const string fromProcString
165  (
166  "procBoundary"
167  + name(procI)
168  + "to"
169  + name(mergedProcI)
170  );
171 
172  if (pp.name() == fromProcString)
173  {
174  isConnected = true;
175  break;
176  }
177  }
178 
179  if (isConnected)
180  {
181  label meshFaceI = pp.start();
182  forAll(pp, i)
183  {
184  addFaces.append(meshFaceI++);
185  }
186  }
187  }
188  }
189  addFaces.shrink();
190 
192  (
193  new faceCoupleInfo
194  (
195  masterMesh,
196  masterFaces,
197  meshToAdd,
198  addFaces,
199  mergeDist, // Absolute merging distance
200  true, // Matching faces identical?
201  false, // If perfect match are faces already ordered
202  // (e.g. processor patches)
203  false // are faces each on separate patch?
204  )
205  );
206  }
207 }
208 
209 
211 (
212  const scalar mergeDist,
213  polyMesh& mesh,
214  labelListList& pointProcAddressing
215 )
216 {
217  // Find out which sets of points get merged and create a map from
218  // mesh point to unique point.
219  Map<label> pointToMaster
220  (
222  (
223  mesh,
224  mergeDist
225  )
226  );
227 
228  Info<< "mergeSharedPoints : detected " << pointToMaster.size()
229  << " points that are to be merged." << endl;
230 
231  if (returnReduce(pointToMaster.size(), sumOp<label>()) == 0)
232  {
233  return autoPtr<mapPolyMesh>(NULL);
234  }
235 
236  polyTopoChange meshMod(mesh);
237 
238  fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
239 
240  // Change the mesh (no inflation). Note: parallel comms allowed.
241  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
242 
243  // Update fields. No inflation, parallel sync.
244  mesh.updateMesh(map);
245 
246  // pointProcAddressing give indices into the master mesh so adapt them
247  // for changed point numbering.
248 
249  // Adapt constructMaps for merged points.
250  forAll(pointProcAddressing, procI)
251  {
252  labelList& constructMap = pointProcAddressing[procI];
253 
254  forAll(constructMap, i)
255  {
256  label oldPointI = constructMap[i];
257 
258  // New label of point after changeMesh.
259  label newPointI = map().reversePointMap()[oldPointI];
260 
261  if (newPointI < -1)
262  {
263  constructMap[i] = -newPointI-2;
264  }
265  else if (newPointI >= 0)
266  {
267  constructMap[i] = newPointI;
268  }
269  else
270  {
272  << "Problem. oldPointI:" << oldPointI
273  << " newPointI:" << newPointI << abort(FatalError);
274  }
275  }
276  }
277 
278  return map;
279 }
280 
281 
283 (
284  const argList& args,
285  const PtrList<Time>& databases,
286  const word& regionDir
287 )
288 {
290 
291  forAll(databases, procI)
292  {
293  fileName pointsInstance
294  (
295  databases[procI].findInstance
296  (
297  regionDir/polyMesh::meshSubDir,
298  "points"
299  )
300  );
301 
302  if (pointsInstance != databases[procI].timeName())
303  {
305  << "Your time was specified as " << databases[procI].timeName()
306  << " but there is no polyMesh/points in that time." << endl
307  << "(there is a points file in " << pointsInstance
308  << ")" << endl
309  << "Please rerun with the correct time specified"
310  << " (through the -constant, -time or -latestTime "
311  << "(at your option)."
312  << endl << exit(FatalError);
313  }
314 
315  Info<< "Reading points from "
316  << databases[procI].caseName()
317  << " for time = " << databases[procI].timeName()
318  << nl << endl;
319 
321  (
322  IOobject
323  (
324  "points",
325  databases[procI].findInstance
326  (
327  regionDir/polyMesh::meshSubDir,
328  "points"
329  ),
330  regionDir/polyMesh::meshSubDir,
331  databases[procI],
334  false
335  )
336  );
337 
338  boundBox domainBb(points, false);
339 
340  bb.min() = min(bb.min(), domainBb.min());
341  bb.max() = max(bb.max(), domainBb.max());
342  }
343 
344  return bb;
345 }
346 
347 
349 (
350  Time& runTime,
351  const fvMesh& masterMesh,
352  const labelListList& cellProcAddressing
353 
354 )
355 {
356  // Write the decomposition as labelList for use with 'manual'
357  // decomposition method.
358  labelIOList cellDecomposition
359  (
360  IOobject
361  (
362  "cellDecomposition",
363  masterMesh.facesInstance(),
364  masterMesh,
367  false
368  ),
369  masterMesh.nCells()
370  );
371 
372  forAll(cellProcAddressing, procI)
373  {
374  const labelList& pCells = cellProcAddressing[procI];
375  UIndirectList<label>(cellDecomposition, pCells) = procI;
376  }
377 
378  cellDecomposition.write();
379 
380  Info<< nl << "Wrote decomposition to "
381  << cellDecomposition.objectPath()
382  << " for use in manual decomposition." << endl;
383 
384 
385  // Write as volScalarField for postprocessing. Change time to 0
386  // if was 'constant'
387  {
388  const scalar oldTime = runTime.value();
389  const label oldIndex = runTime.timeIndex();
390  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
391  {
392  runTime.setTime(0, oldIndex+1);
393  }
394 
395  volScalarField cellDist
396  (
397  IOobject
398  (
399  "cellDist",
400  runTime.timeName(),
401  masterMesh,
404  ),
405  masterMesh,
406  dimensionedScalar("cellDist", dimless, 0),
407  zeroGradientFvPatchScalarField::typeName
408  );
409 
410  forAll(cellDecomposition, cellI)
411  {
412  cellDist[cellI] = cellDecomposition[cellI];
413  }
414 
415  cellDist.write();
416 
417  Info<< nl << "Wrote decomposition as volScalarField to "
418  << cellDist.name() << " for use in postprocessing."
419  << endl;
420 
421  // Restore time
422  runTime.setTime(oldTime, oldIndex);
423  }
424 }
425 
426 
427 int main(int argc, char *argv[])
428 {
430  (
431  "reconstruct a mesh using geometric information only"
432  );
433 
434  // Enable -constant ... if someone really wants it
435  // Enable -withZero to prevent accidentally trashing the initial fields
436  timeSelector::addOptions(true, true);
439  (
440  "mergeTol",
441  "scalar",
442  "specify the merge distance relative to the bounding box size "
443  "(default 1e-7)"
444  );
446  (
447  "fullMatch",
448  "do (slower) geometric matching on all boundary faces"
449  );
451  (
452  "cellDist",
453  "write cell distribution as a labelList - for use with 'manual' "
454  "decomposition method or as a volScalarField for post-processing."
455  );
456 
457  #include "addRegionOption.H"
458  #include "setRootCase.H"
459  #include "createTime.H"
460 
461  Info<< "This is an experimental tool which tries to merge"
462  << " individual processor" << nl
463  << "meshes back into one master mesh. Use it if the original"
464  << " master mesh has" << nl
465  << "been deleted or if the processor meshes have been modified"
466  << " (topology change)." << nl
467  << "This tool will write the resulting mesh to a new time step"
468  << " and construct" << nl
469  << "xxxxProcAddressing files in the processor meshes so"
470  << " reconstructPar can be" << nl
471  << "used to regenerate the fields on the master mesh." << nl
472  << nl
473  << "Not well tested & use at your own risk!" << nl
474  << endl;
475 
476 
478  word regionDir = word::null;
479 
480  if
481  (
484  )
485  {
486  regionDir = regionName;
487  Info<< "Operating on region " << regionName << nl << endl;
488  }
489 
490  scalar mergeTol = defaultMergeTol;
491  args.optionReadIfPresent("mergeTol", mergeTol);
492 
493  scalar writeTol = Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
494 
495  Info<< "Merge tolerance : " << mergeTol << nl
496  << "Write tolerance : " << writeTol << endl;
497 
498  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
499  {
501  << "Your current settings specify ASCII writing with "
502  << IOstream::defaultPrecision() << " digits precision." << endl
503  << "Your merging tolerance (" << mergeTol << ") is finer than this."
504  << endl
505  << "Please change your writeFormat to binary"
506  << " or increase the writePrecision" << endl
507  << "or adjust the merge tolerance (-mergeTol)."
508  << exit(FatalError);
509  }
510 
511 
512  const bool fullMatch = args.optionFound("fullMatch");
513 
514  if (fullMatch)
515  {
516  Info<< "Doing geometric matching on all boundary faces." << nl << endl;
517  }
518  else
519  {
520  Info<< "Doing geometric matching on correct procBoundaries only."
521  << nl << "This assumes a correct decomposition." << endl;
522  }
523 
524  bool writeCellDist = args.optionFound("cellDist");
525 
526 
527  int nProcs = 0;
528 
529  while
530  (
531  isDir
532  (
533  args.rootPath()
534  / args.caseName()
535  / fileName(word("processor") + name(nProcs))
536  )
537  )
538  {
539  nProcs++;
540  }
541 
542  Info<< "Found " << nProcs << " processor directories" << nl << endl;
543 
544 
545  // Read all time databases
546  PtrList<Time> databases(nProcs);
547 
548  forAll(databases, procI)
549  {
550  Info<< "Reading database "
551  << args.caseName()/fileName(word("processor") + name(procI))
552  << endl;
553 
554  databases.set
555  (
556  procI,
557  new Time
558  (
560  args.rootPath(),
561  args.caseName()/fileName(word("processor") + name(procI))
562  )
563  );
564  }
565 
566  // Use the times list from the master processor
567  // and select a subset based on the command-line options
569  (
570  databases[0].times(),
571  args
572  );
573 
574  // Loop over all times
575  forAll(timeDirs, timeI)
576  {
577  // Set time for global database
578  runTime.setTime(timeDirs[timeI], timeI);
579 
580  Info<< "Time = " << runTime.timeName() << nl << endl;
581 
582  // Set time for all databases
583  forAll(databases, procI)
584  {
585  databases[procI].setTime(timeDirs[timeI], timeI);
586  }
587 
588  const fileName meshPath =
589  databases[0].path()
590  /databases[0].timeName()
591  /regionDir
593 
594  if (!isFile(meshPath/"faces"))
595  {
596  Info<< "No mesh." << nl << endl;
597  continue;
598  }
599 
600 
601  // Read point on individual processors to determine merge tolerance
602  // (otherwise single cell domains might give problems)
603 
604  const boundBox bb = procBounds(args, databases, regionDir);
605  const scalar mergeDist = mergeTol*bb.mag();
606 
607  Info<< "Overall mesh bounding box : " << bb << nl
608  << "Relative tolerance : " << mergeTol << nl
609  << "Absolute matching distance : " << mergeDist << nl
610  << endl;
611 
612 
613  // Addressing from processor to reconstructed case
614  labelListList cellProcAddressing(nProcs);
616  labelListList pointProcAddressing(nProcs);
617  labelListList boundaryProcAddressing(nProcs);
618 
619  // Internal faces on the final reconstructed mesh
620  label masterInternalFaces;
621 
622  // Owner addressing on the final reconstructed mesh
623  labelList masterOwner;
624 
625  {
626  // Construct empty mesh.
627  Info<< "Constructing empty mesh to add to." << nl << endl;
628  fvMesh masterMesh
629  (
630  IOobject
631  (
632  regionName,
633  runTime.timeName(),
634  runTime,
636  ),
637  xferCopy(pointField()),
638  xferCopy(faceList()),
639  xferCopy(cellList())
640  );
641 
642  for (label procI = 0; procI < nProcs; procI++)
643  {
644  Info<< "Reading mesh to add from "
645  << databases[procI].caseName()
646  << " for time = " << databases[procI].timeName()
647  << nl << endl;
648 
649  fvMesh meshToAdd
650  (
651  IOobject
652  (
653  regionName,
654  databases[procI].timeName(),
655  databases[procI]
656  )
657  );
658 
659  // Initialize its addressing
660  cellProcAddressing[procI] = identity(meshToAdd.nCells());
661  faceProcAddressing[procI] = identity(meshToAdd.nFaces());
662  pointProcAddressing[procI] = identity(meshToAdd.nPoints());
663  boundaryProcAddressing[procI] =
664  identity(meshToAdd.boundaryMesh().size());
665 
666 
667  // Find geometrically shared points/faces.
669  (
670  fullMatch,
671  procI,
672  masterMesh,
673  meshToAdd,
674  mergeDist
675  );
676 
677 
678  // Add elements to mesh
679  Info<< "Adding to master mesh" << nl << endl;
680 
682  (
683  masterMesh,
684  meshToAdd,
685  couples
686  );
687 
688  // Update all addressing so xxProcAddressing points to correct
689  // item in masterMesh.
690 
691  // Processors that were already in masterMesh
692  for (label mergedI = 0; mergedI < procI; mergedI++)
693  {
694  renumber(map().oldCellMap(), cellProcAddressing[mergedI]);
695  renumber(map().oldFaceMap(), faceProcAddressing[mergedI]);
696  renumber(map().oldPointMap(), pointProcAddressing[mergedI]);
697  // Note: boundary is special since can contain -1.
698  renumber
699  (
700  map().oldPatchMap(),
701  boundaryProcAddressing[mergedI]
702  );
703  }
704 
705  // Added processor
706  renumber(map().addedCellMap(), cellProcAddressing[procI]);
707  renumber(map().addedFaceMap(), faceProcAddressing[procI]);
708  renumber(map().addedPointMap(), pointProcAddressing[procI]);
709  renumber(map().addedPatchMap(), boundaryProcAddressing[procI]);
710 
711  Info<< endl;
712  }
713 
714  // See if any points on the mastermesh have become connected
715  // because of connections through processor meshes.
716  mergeSharedPoints(mergeDist, masterMesh, pointProcAddressing);
717 
718  // Save some properties on the reconstructed mesh
719  masterInternalFaces = masterMesh.nInternalFaces();
720  masterOwner = masterMesh.faceOwner();
721 
722 
723  Info<< "\nWriting merged mesh to "
724  << runTime.path()/runTime.timeName()
725  << nl << endl;
726 
727  if (!masterMesh.write())
728  {
730  << "Failed writing polyMesh."
731  << exit(FatalError);
732  }
733 
734  if (writeCellDist)
735  {
736  writeCellDistance(runTime, masterMesh, cellProcAddressing);
737  }
738  }
739 
740 
741  // Write the addressing
742 
743  Info<< "Reconstructing the addressing from the processor meshes"
744  << " to the newly reconstructed mesh" << nl << endl;
745 
746  forAll(databases, procI)
747  {
748  Info<< "Reading processor " << procI << " mesh from "
749  << databases[procI].caseName() << endl;
750 
751  polyMesh procMesh
752  (
753  IOobject
754  (
755  regionName,
756  databases[procI].timeName(),
757  databases[procI]
758  )
759  );
760 
761 
762  // From processor point to reconstructed mesh point
763 
764  Info<< "Writing pointProcAddressing to "
765  << databases[procI].caseName()
766  /procMesh.facesInstance()
768  << endl;
769 
771  (
772  IOobject
773  (
774  "pointProcAddressing",
775  procMesh.facesInstance(),
777  procMesh,
780  false // Do not register
781  ),
782  pointProcAddressing[procI]
783  ).write();
784 
785 
786  // From processor face to reconstructed mesh face
787 
788  Info<< "Writing faceProcAddressing to "
789  << databases[procI].caseName()
790  /procMesh.facesInstance()
792  << endl;
793 
794  labelIOList faceProcAddr
795  (
796  IOobject
797  (
798  "faceProcAddressing",
799  procMesh.facesInstance(),
801  procMesh,
804  false // Do not register
805  ),
806  faceProcAddressing[procI]
807  );
808 
809  // Now add turning index to faceProcAddressing.
810  // See reconstructPar for meaning of turning index.
811  forAll(faceProcAddr, procFaceI)
812  {
813  label masterFaceI = faceProcAddr[procFaceI];
814 
815  if
816  (
817  !procMesh.isInternalFace(procFaceI)
818  && masterFaceI < masterInternalFaces
819  )
820  {
821  // proc face is now external but used to be internal face.
822  // Check if we have owner or neighbour.
823 
824  label procOwn = procMesh.faceOwner()[procFaceI];
825  label masterOwn = masterOwner[masterFaceI];
826 
827  if (cellProcAddressing[procI][procOwn] == masterOwn)
828  {
829  // No turning. Offset by 1.
830  faceProcAddr[procFaceI]++;
831  }
832  else
833  {
834  // Turned face.
835  faceProcAddr[procFaceI] =
836  -1 - faceProcAddr[procFaceI];
837  }
838  }
839  else
840  {
841  // No turning. Offset by 1.
842  faceProcAddr[procFaceI]++;
843  }
844  }
845 
846  faceProcAddr.write();
847 
848 
849  // From processor cell to reconstructed mesh cell
850 
851  Info<< "Writing cellProcAddressing to "
852  << databases[procI].caseName()
853  /procMesh.facesInstance()
855  << endl;
856 
858  (
859  IOobject
860  (
861  "cellProcAddressing",
862  procMesh.facesInstance(),
864  procMesh,
867  false // Do not register
868  ),
869  cellProcAddressing[procI]
870  ).write();
871 
872 
873 
874  // From processor patch to reconstructed mesh patch
875 
876  Info<< "Writing boundaryProcAddressing to "
877  << databases[procI].caseName()
878  /procMesh.facesInstance()
880  << endl;
881 
883  (
884  IOobject
885  (
886  "boundaryProcAddressing",
887  procMesh.facesInstance(),
889  procMesh,
892  false // Do not register
893  ),
894  boundaryProcAddressing[procI]
895  ).write();
896 
897  Info<< endl;
898  }
899  }
900 
901 
902  Info<< "\nEnd\n" << endl;
903 
904  return 0;
905 }
906 
907 
908 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Vector::min
static const Vector min
Definition: Vector.H:83
Foam::boundBox::mag
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField< vector >
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< label >
mergeSharedPoints
autoPtr< mapPolyMesh > mergeSharedPoints(const scalar mergeDist, polyMesh &mesh, labelListList &pointProcAddressing)
Definition: reconstructParMesh.C:208
Foam::Time::writeFormat
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:303
Foam::fvMeshAdder::add
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:74
Foam::boundBox::invertedBox
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- VGREAT.
Definition: boundBox.H:79
Foam::Vector::max
static const Vector max
Definition: Vector.H:82
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
polyTopoChange.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::Map< label >
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:97
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
oldTime
rho oldTime()
Foam::argList::rootPath
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
fvMeshAdder.H
Foam::renumber
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:32
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
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:66
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
addRegionOption.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::cellList
List< cell > cellList
list of cells
Definition: cellList.H:42
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
mapAddedPolyMesh.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
faceCoupleInfo.H
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
determineCoupledFaces
autoPtr< faceCoupleInfo > determineCoupledFaces(const bool fullMatch, const label procI, const polyMesh &masterMesh, const polyMesh &meshToAdd, const scalar mergeDist)
Definition: reconstructParMesh.C:81
Foam::FatalError
error FatalError
processorPolyPatch.H
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::Time::setTime
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:964
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
procBounds
boundBox procBounds(const argList &args, const PtrList< Time > &databases, const word &regionDir)
Definition: reconstructParMesh.C:280
Foam::isFile
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:622
Foam::isDir
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:615
polyMeshAdder.H
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::polyMeshAdder::findSharedPoints
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
Definition: polyMeshAdder.C:2000
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::sumOp
Definition: ops.H:162
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:281
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Foam::xferCopy
Xfer< T > xferCopy(const T &)
Construct by copying the contents of the arg.
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
labelIOList.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::IOList< label >
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
timeSelector.H
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
faceProcAddressing
PtrList< labelIOList > & faceProcAddressing
Definition: checkFaceAddressingComp.H:9
Foam::faceCoupleInfo
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
Definition: faceCoupleInfo.H:157
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::timeSelector::select
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::argList::caseName
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:161
defaultMergeTol
static const scalar defaultMergeTol
Definition: reconstructParMesh.C:57
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
Foam::Time::controlDictName
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:213
zeroGradientFvPatchFields.H
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
main
int main(int argc, char *argv[])
Definition: reconstructParMesh.C:424
Foam::polyMeshAdder::mergePoints
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
Definition: polyMeshAdder.C:2212
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
writeCellDistance
void writeCellDistance(Time &runTime, const fvMesh &masterMesh, const labelListList &cellProcAddressing)
Definition: reconstructParMesh.C:346