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-2016 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 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 word fromProcString
165  (
166  processorPolyPatch::newName(proci, mergedProci)
167  );
168 
169  if (pp.name() == fromProcString)
170  {
171  isConnected = true;
172  break;
173  }
174  }
175 
176  if (isConnected)
177  {
178  label meshFacei = pp.start();
179  forAll(pp, i)
180  {
181  addFaces.append(meshFacei++);
182  }
183  }
184  }
185  }
186  addFaces.shrink();
187 
189  (
190  new faceCoupleInfo
191  (
192  masterMesh,
193  masterFaces,
194  meshToAdd,
195  addFaces,
196  mergeDist, // Absolute merging distance
197  true, // Matching faces identical?
198  false, // If perfect match are faces already ordered
199  // (e.g. processor patches)
200  false // are faces each on separate patch?
201  )
202  );
203  }
204 }
205 
206 
207 autoPtr<mapPolyMesh> mergeSharedPoints
208 (
209  const scalar mergeDist,
210  polyMesh& mesh,
211  labelListList& pointProcAddressing
212 )
213 {
214  // Find out which sets of points get merged and create a map from
215  // mesh point to unique point.
216  Map<label> pointToMaster
217  (
219  (
220  mesh,
221  mergeDist
222  )
223  );
224 
225  Info<< "mergeSharedPoints : detected " << pointToMaster.size()
226  << " points that are to be merged." << endl;
227 
228  if (returnReduce(pointToMaster.size(), sumOp<label>()) == 0)
229  {
230  return autoPtr<mapPolyMesh>(NULL);
231  }
232 
233  polyTopoChange meshMod(mesh);
234 
235  fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
236 
237  // Change the mesh (no inflation). Note: parallel comms allowed.
238  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
239 
240  // Update fields. No inflation, parallel sync.
241  mesh.updateMesh(map);
242 
243  // pointProcAddressing give indices into the master mesh so adapt them
244  // for changed point numbering.
245 
246  // Adapt constructMaps for merged points.
247  forAll(pointProcAddressing, proci)
248  {
249  labelList& constructMap = pointProcAddressing[proci];
250 
251  forAll(constructMap, i)
252  {
253  label oldPointi = constructMap[i];
254 
255  // New label of point after changeMesh.
256  label newPointi = map().reversePointMap()[oldPointi];
257 
258  if (newPointi < -1)
259  {
260  constructMap[i] = -newPointi-2;
261  }
262  else if (newPointi >= 0)
263  {
264  constructMap[i] = newPointi;
265  }
266  else
267  {
269  << "Problem. oldPointi:" << oldPointi
270  << " newPointi:" << newPointi << abort(FatalError);
271  }
272  }
273  }
274 
275  return map;
276 }
277 
278 
279 boundBox procBounds
280 (
281  const argList& args,
282  const PtrList<Time>& databases,
283  const word& regionDir
284 )
285 {
287 
288  forAll(databases, proci)
289  {
290  fileName pointsInstance
291  (
292  databases[proci].findInstance
293  (
294  regionDir/polyMesh::meshSubDir,
295  "points"
296  )
297  );
298 
299  if (pointsInstance != databases[proci].timeName())
300  {
302  << "Your time was specified as " << databases[proci].timeName()
303  << " but there is no polyMesh/points in that time." << endl
304  << "(there is a points file in " << pointsInstance
305  << ")" << endl
306  << "Please rerun with the correct time specified"
307  << " (through the -constant, -time or -latestTime "
308  << "(at your option)."
309  << endl << exit(FatalError);
310  }
311 
312  Info<< "Reading points from "
313  << databases[proci].caseName()
314  << " for time = " << databases[proci].timeName()
315  << nl << endl;
316 
318  (
319  IOobject
320  (
321  "points",
322  databases[proci].findInstance
323  (
324  regionDir/polyMesh::meshSubDir,
325  "points"
326  ),
327  regionDir/polyMesh::meshSubDir,
328  databases[proci],
330  IOobject::NO_WRITE,
331  false
332  )
333  );
334 
335  boundBox domainBb(points, false);
336 
337  bb.min() = min(bb.min(), domainBb.min());
338  bb.max() = max(bb.max(), domainBb.max());
339  }
340 
341  return bb;
342 }
343 
344 
345 void writeCellDistance
346 (
347  Time& runTime,
348  const fvMesh& masterMesh,
349  const labelListList& cellProcAddressing
350 
351 )
352 {
353  // Write the decomposition as labelList for use with 'manual'
354  // decomposition method.
355  labelIOList cellDecomposition
356  (
357  IOobject
358  (
359  "cellDecomposition",
360  masterMesh.facesInstance(),
361  masterMesh,
364  false
365  ),
366  masterMesh.nCells()
367  );
368 
369  forAll(cellProcAddressing, proci)
370  {
371  const labelList& pCells = cellProcAddressing[proci];
372  UIndirectList<label>(cellDecomposition, pCells) = proci;
373  }
374 
375  cellDecomposition.write();
376 
377  Info<< nl << "Wrote decomposition to "
378  << cellDecomposition.objectPath()
379  << " for use in manual decomposition." << endl;
380 
381 
382  // Write as volScalarField for postprocessing. Change time to 0
383  // if was 'constant'
384  {
385  const scalar oldTime = runTime.value();
386  const label oldIndex = runTime.timeIndex();
387  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
388  {
389  runTime.setTime(0, oldIndex+1);
390  }
391 
392  volScalarField cellDist
393  (
394  IOobject
395  (
396  "cellDist",
397  runTime.timeName(),
398  masterMesh,
401  ),
402  masterMesh,
403  dimensionedScalar("cellDist", dimless, 0),
404  extrapolatedCalculatedFvPatchScalarField::typeName
405  );
406 
407  forAll(cellDecomposition, celli)
408  {
409  cellDist[celli] = cellDecomposition[celli];
410  }
411  cellDist.correctBoundaryConditions();
412 
413  cellDist.write();
414 
415  Info<< nl << "Wrote decomposition as volScalarField to "
416  << cellDist.name() << " for use in postprocessing."
417  << endl;
418 
419  // Restore time
420  runTime.setTime(oldTime, oldIndex);
421  }
422 }
423 
424 
425 int main(int argc, char *argv[])
426 {
428  (
429  "reconstruct a mesh using geometric information only"
430  );
431 
432  // Enable -constant ... if someone really wants it
433  // Enable -withZero to prevent accidentally trashing the initial fields
434  timeSelector::addOptions(true, true);
437  (
438  "mergeTol",
439  "scalar",
440  "specify the merge distance relative to the bounding box size "
441  "(default 1e-7)"
442  );
444  (
445  "fullMatch",
446  "do (slower) geometric matching on all boundary faces"
447  );
449  (
450  "cellDist",
451  "write cell distribution as a labelList - for use with 'manual' "
452  "decomposition method or as a volScalarField for post-processing."
453  );
454 
455  #include "addRegionOption.H"
456  #include "setRootCase.H"
457  #include "createTime.H"
458 
459  Info<< "This is an experimental tool which tries to merge"
460  << " individual processor" << nl
461  << "meshes back into one master mesh. Use it if the original"
462  << " master mesh has" << nl
463  << "been deleted or if the processor meshes have been modified"
464  << " (topology change)." << nl
465  << "This tool will write the resulting mesh to a new time step"
466  << " and construct" << nl
467  << "xxxxProcAddressing files in the processor meshes so"
468  << " reconstructPar can be" << nl
469  << "used to regenerate the fields on the master mesh." << nl
470  << nl
471  << "Not well tested & use at your own risk!" << nl
472  << endl;
473 
474 
475  word regionName = polyMesh::defaultRegion;
476  word regionDir = word::null;
477 
478  if
479  (
480  args.optionReadIfPresent("region", regionName)
481  && regionName != polyMesh::defaultRegion
482  )
483  {
484  regionDir = regionName;
485  Info<< "Operating on region " << regionName << nl << endl;
486  }
487 
488  scalar mergeTol = defaultMergeTol;
489  args.optionReadIfPresent("mergeTol", mergeTol);
490 
491  scalar writeTol = Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
492 
493  Info<< "Merge tolerance : " << mergeTol << nl
494  << "Write tolerance : " << writeTol << endl;
495 
496  if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
497  {
499  << "Your current settings specify ASCII writing with "
500  << IOstream::defaultPrecision() << " digits precision." << endl
501  << "Your merging tolerance (" << mergeTol << ") is finer than this."
502  << endl
503  << "Please change your writeFormat to binary"
504  << " or increase the writePrecision" << endl
505  << "or adjust the merge tolerance (-mergeTol)."
506  << exit(FatalError);
507  }
508 
509 
510  const bool fullMatch = args.optionFound("fullMatch");
511 
512  if (fullMatch)
513  {
514  Info<< "Doing geometric matching on all boundary faces." << nl << endl;
515  }
516  else
517  {
518  Info<< "Doing geometric matching on correct procBoundaries only."
519  << nl << "This assumes a correct decomposition." << endl;
520  }
521 
522  bool writeCellDist = args.optionFound("cellDist");
523 
524 
525  int nProcs = 0;
526 
527  while
528  (
529  isDir
530  (
531  args.rootPath()
532  / args.caseName()
533  / fileName(word("processor") + name(nProcs))
534  )
535  )
536  {
537  nProcs++;
538  }
539 
540  Info<< "Found " << nProcs << " processor directories" << nl << endl;
541 
542 
543  // Read all time databases
544  PtrList<Time> databases(nProcs);
545 
546  forAll(databases, proci)
547  {
548  Info<< "Reading database "
549  << args.caseName()/fileName(word("processor") + name(proci))
550  << endl;
551 
552  databases.set
553  (
554  proci,
555  new Time
556  (
558  args.rootPath(),
559  args.caseName()/fileName(word("processor") + name(proci))
560  )
561  );
562  }
563 
564  // Use the times list from the master processor
565  // and select a subset based on the command-line options
567  (
568  databases[0].times(),
569  args
570  );
571 
572  // Loop over all times
573  forAll(timeDirs, timeI)
574  {
575  // Set time for global database
576  runTime.setTime(timeDirs[timeI], timeI);
577 
578  Info<< "Time = " << runTime.timeName() << nl << endl;
579 
580  // Set time for all databases
581  forAll(databases, proci)
582  {
583  databases[proci].setTime(timeDirs[timeI], timeI);
584  }
585 
586  const fileName meshPath =
587  databases[0].path()
588  /databases[0].timeName()
589  /regionDir
591 
592  if (!isFile(meshPath/"faces"))
593  {
594  Info<< "No mesh." << nl << endl;
595  continue;
596  }
597 
598 
599  // Read point on individual processors to determine merge tolerance
600  // (otherwise single cell domains might give problems)
601 
602  const boundBox bb = procBounds(args, databases, regionDir);
603  const scalar mergeDist = mergeTol*bb.mag();
604 
605  Info<< "Overall mesh bounding box : " << bb << nl
606  << "Relative tolerance : " << mergeTol << nl
607  << "Absolute matching distance : " << mergeDist << nl
608  << endl;
609 
610 
611  // Addressing from processor to reconstructed case
612  labelListList cellProcAddressing(nProcs);
614  labelListList pointProcAddressing(nProcs);
615  labelListList boundaryProcAddressing(nProcs);
616 
617  // Internal faces on the final reconstructed mesh
618  label masterInternalFaces;
619 
620  // Owner addressing on the final reconstructed mesh
621  labelList masterOwner;
622 
623  {
624  // Construct empty mesh.
625  Info<< "Constructing empty mesh to add to." << nl << endl;
626  fvMesh masterMesh
627  (
628  IOobject
629  (
630  regionName,
631  runTime.timeName(),
632  runTime,
634  ),
635  xferCopy(pointField()),
636  xferCopy(faceList()),
637  xferCopy(cellList())
638  );
639 
640  for (label proci = 0; proci < nProcs; proci++)
641  {
642  Info<< "Reading mesh to add from "
643  << databases[proci].caseName()
644  << " for time = " << databases[proci].timeName()
645  << nl << endl;
646 
647  fvMesh meshToAdd
648  (
649  IOobject
650  (
651  regionName,
652  databases[proci].timeName(),
653  databases[proci]
654  )
655  );
656 
657  // Initialize its addressing
658  cellProcAddressing[proci] = identity(meshToAdd.nCells());
659  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
660  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
661  boundaryProcAddressing[proci] =
662  identity(meshToAdd.boundaryMesh().size());
663 
664 
665  // Find geometrically shared points/faces.
666  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
667  (
668  fullMatch,
669  proci,
670  masterMesh,
671  meshToAdd,
672  mergeDist
673  );
674 
675 
676  // Add elements to mesh
677  Info<< "Adding to master mesh" << nl << endl;
678 
680  (
681  masterMesh,
682  meshToAdd,
683  couples
684  );
685 
686  // Update all addressing so xxProcAddressing points to correct
687  // item in masterMesh.
688 
689  // Processors that were already in masterMesh
690  for (label mergedI = 0; mergedI < proci; mergedI++)
691  {
692  renumber(map().oldCellMap(), cellProcAddressing[mergedI]);
693  renumber(map().oldFaceMap(), faceProcAddressing[mergedI]);
694  renumber(map().oldPointMap(), pointProcAddressing[mergedI]);
695  // Note: boundary is special since can contain -1.
696  renumber
697  (
698  map().oldPatchMap(),
699  boundaryProcAddressing[mergedI]
700  );
701  }
702 
703  // Added processor
704  renumber(map().addedCellMap(), cellProcAddressing[proci]);
705  renumber(map().addedFaceMap(), faceProcAddressing[proci]);
706  renumber(map().addedPointMap(), pointProcAddressing[proci]);
707  renumber(map().addedPatchMap(), boundaryProcAddressing[proci]);
708 
709  Info<< endl;
710  }
711 
712  // See if any points on the mastermesh have become connected
713  // because of connections through processor meshes.
714  mergeSharedPoints(mergeDist, masterMesh, pointProcAddressing);
715 
716  // Save some properties on the reconstructed mesh
717  masterInternalFaces = masterMesh.nInternalFaces();
718  masterOwner = masterMesh.faceOwner();
719 
720 
721  Info<< "\nWriting merged mesh to "
722  << runTime.path()/runTime.timeName()
723  << nl << endl;
724 
725  if (!masterMesh.write())
726  {
728  << "Failed writing polyMesh."
729  << exit(FatalError);
730  }
731 
732  if (writeCellDist)
733  {
734  writeCellDistance(runTime, masterMesh, cellProcAddressing);
735  }
736  }
737 
738 
739  // Write the addressing
740 
741  Info<< "Reconstructing the addressing from the processor meshes"
742  << " to the newly reconstructed mesh" << nl << endl;
743 
744  forAll(databases, proci)
745  {
746  Info<< "Reading processor " << proci << " mesh from "
747  << databases[proci].caseName() << endl;
748 
749  polyMesh procMesh
750  (
751  IOobject
752  (
753  regionName,
754  databases[proci].timeName(),
755  databases[proci]
756  )
757  );
758 
759 
760  // From processor point to reconstructed mesh point
761 
762  Info<< "Writing pointProcAddressing to "
763  << databases[proci].caseName()
764  /procMesh.facesInstance()
766  << endl;
767 
769  (
770  IOobject
771  (
772  "pointProcAddressing",
773  procMesh.facesInstance(),
775  procMesh,
778  false // Do not register
779  ),
780  pointProcAddressing[proci]
781  ).write();
782 
783 
784  // From processor face to reconstructed mesh face
785 
786  Info<< "Writing faceProcAddressing to "
787  << databases[proci].caseName()
788  /procMesh.facesInstance()
790  << endl;
791 
792  labelIOList faceProcAddr
793  (
794  IOobject
795  (
796  "faceProcAddressing",
797  procMesh.facesInstance(),
799  procMesh,
802  false // Do not register
803  ),
804  faceProcAddressing[proci]
805  );
806 
807  // Now add turning index to faceProcAddressing.
808  // See reconstructPar for meaning of turning index.
809  forAll(faceProcAddr, procFacei)
810  {
811  label masterFacei = faceProcAddr[procFacei];
812 
813  if
814  (
815  !procMesh.isInternalFace(procFacei)
816  && masterFacei < masterInternalFaces
817  )
818  {
819  // proc face is now external but used to be internal face.
820  // Check if we have owner or neighbour.
821 
822  label procOwn = procMesh.faceOwner()[procFacei];
823  label masterOwn = masterOwner[masterFacei];
824 
825  if (cellProcAddressing[proci][procOwn] == masterOwn)
826  {
827  // No turning. Offset by 1.
828  faceProcAddr[procFacei]++;
829  }
830  else
831  {
832  // Turned face.
833  faceProcAddr[procFacei] =
834  -1 - faceProcAddr[procFacei];
835  }
836  }
837  else
838  {
839  // No turning. Offset by 1.
840  faceProcAddr[procFacei]++;
841  }
842  }
843 
844  faceProcAddr.write();
845 
846 
847  // From processor cell to reconstructed mesh cell
848 
849  Info<< "Writing cellProcAddressing to "
850  << databases[proci].caseName()
851  /procMesh.facesInstance()
853  << endl;
854 
856  (
857  IOobject
858  (
859  "cellProcAddressing",
860  procMesh.facesInstance(),
862  procMesh,
865  false // Do not register
866  ),
867  cellProcAddressing[proci]
868  ).write();
869 
870 
871 
872  // From processor patch to reconstructed mesh patch
873 
874  Info<< "Writing boundaryProcAddressing to "
875  << databases[proci].caseName()
876  /procMesh.facesInstance()
878  << endl;
879 
881  (
882  IOobject
883  (
884  "boundaryProcAddressing",
885  procMesh.facesInstance(),
887  procMesh,
890  false // Do not register
891  ),
892  boundaryProcAddressing[proci]
893  ).write();
894 
895  Info<< endl;
896  }
897  }
898 
899 
900  Info<< "End.\n" << endl;
901 
902  return 0;
903 }
904 
905 
906 // ************************************************************************* //
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
static const Form max
Definition: VectorSpace.H:112
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const double e
Elementary charge.
Definition: doubleFloat.H:78
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
error FatalError
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
PtrList< labelIOList > & faceProcAddressing
Foam::word regionName
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
List< face > faceList
Definition: faceListFwd.H:43
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:492
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
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
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
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.
const Type & value() const
Return const reference to value.
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
bool isDir(const fileName &)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:486
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
static const Form min
Definition: VectorSpace.H:113
const fileName & rootPath() const
Return root path.
Definition: argListI.H:36
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
label nCells() const
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:42
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
const pointField & points
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:93
static const word null
An empty word.
Definition: word.H:77
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
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:924
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:291
const word & name() const
Return name.
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.
static const char nl
Definition: Ostream.H:262
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label nFaces() const
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- VGREAT.
Definition: boundBox.H:82
label patchi
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:62
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.
virtual bool write() const
Write using setting from DB.
label nPoints() const
messageStream Info
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:53
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:83
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
fileName path() const
Return path.
Definition: Time.H:269
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
label nInternalFaces() const
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
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.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29