polyMeshFromShapeMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 #include "Time.H"
28 #include "primitiveMesh.H"
29 #include "DynamicList.H"
30 #include "indexedOctree.H"
31 #include "treeDataCell.H"
32 #include "globalMeshData.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 Foam::labelListList Foam::polyMesh::cellShapePointCells
37 (
38  const cellShapeList& c
39 ) const
40 {
41  List<DynamicList<label, primitiveMesh::cellsPerPoint_>>
42  pc(points().size());
43 
44  // For each cell
45  forAll(c, i)
46  {
47  // For each vertex
48  const labelList& labels = c[i];
49 
50  forAll(labels, j)
51  {
52  // Set working point label
53  label curPoint = labels[j];
54  DynamicList<label, primitiveMesh::cellsPerPoint_>& curPointCells =
55  pc[curPoint];
56 
57  // Enter the cell label in the point's cell list
58  curPointCells.append(i);
59  }
60  }
61 
62  labelListList pointCellAddr(pc.size());
63 
64  forAll(pc, pointi)
65  {
66  pointCellAddr[pointi].transfer(pc[pointi]);
67  }
68 
69  return pointCellAddr;
70 }
71 
72 
73 Foam::labelList Foam::polyMesh::facePatchFaceCells
74 (
75  const faceList& patchFaces,
76  const labelListList& pointCells,
77  const faceListList& cellsFaceShapes,
78  const label patchID
79 ) const
80 {
81  bool found;
82 
83  labelList FaceCells(patchFaces.size());
84 
85  forAll(patchFaces, fI)
86  {
87  found = false;
88 
89  const face& curFace = patchFaces[fI];
90  const labelList& facePoints = patchFaces[fI];
91 
92  forAll(facePoints, pointi)
93  {
94  const labelList& facePointCells = pointCells[facePoints[pointi]];
95 
96  forAll(facePointCells, celli)
97  {
98  faceList cellFaces = cellsFaceShapes[facePointCells[celli]];
99 
100  forAll(cellFaces, cellFace)
101  {
102  if (face::sameVertices(cellFaces[cellFace], curFace))
103  {
104  // Found the cell corresponding to this face
105  FaceCells[fI] = facePointCells[celli];
106 
107  found = true;
108  }
109  if (found) break;
110  }
111  if (found) break;
112  }
113  if (found) break;
114  }
115 
116  if (!found)
117  {
119  << "face " << fI << " in patch " << patchID
120  << " does not have neighbour cell"
121  << " face: " << patchFaces[fI]
122  << abort(FatalError);
123  }
124  }
125 
126  return FaceCells;
127 }
128 
129 
130 void Foam::polyMesh::setTopology
131 (
132  const cellShapeList& cellsAsShapes,
133  const faceListList& boundaryFaces,
134  const wordList& boundaryPatchNames,
135  labelList& patchSizes,
136  labelList& patchStarts,
137  label& defaultPatchStart,
138  label& nFaces,
139  cellList& cells
140 )
141 {
142  // Calculate the faces of all cells
143  // Initialise maximum possible number of mesh faces to 0
144  label maxFaces = 0;
145 
146  // Set up a list of face shapes for each cell
147  faceListList cellsFaceShapes(cellsAsShapes.size());
148  cells.setSize(cellsAsShapes.size());
149 
150  forAll(cellsFaceShapes, celli)
151  {
152  cellsFaceShapes[celli] = cellsAsShapes[celli].faces();
153 
154  cells[celli].setSize(cellsFaceShapes[celli].size());
155 
156  // Initialise cells to -1 to flag undefined faces
157  static_cast<labelList&>(cells[celli]) = -1;
158 
159  // Count maximum possible number of mesh faces
160  maxFaces += cellsFaceShapes[celli].size();
161  }
162 
163  // Set size of faces array to maximum possible number of mesh faces
164  faces_.setSize(maxFaces);
165 
166  // Initialise number of faces to 0
167  nFaces = 0;
168 
169  // Set reference to point-cell addressing
170  labelListList PointCells = cellShapePointCells(cellsAsShapes);
171 
172  bool found = false;
173 
174  forAll(cells, celli)
175  {
176  // Note:
177  // Insertion cannot be done in one go as the faces need to be
178  // added into the list in the increasing order of neighbour
179  // cells. Therefore, all neighbours will be detected first
180  // and then added in the correct order.
181 
182  const faceList& curFaces = cellsFaceShapes[celli];
183 
184  // Record the neighbour cell
185  labelList neiCells(curFaces.size(), -1);
186 
187  // Record the face of neighbour cell
188  labelList faceOfNeiCell(curFaces.size(), -1);
189 
190  label nNeighbours = 0;
191 
192  // For all faces ...
193  forAll(curFaces, facei)
194  {
195  // Skip faces that have already been matched
196  if (cells[celli][facei] >= 0) continue;
197 
198  found = false;
199 
200  const face& curFace = curFaces[facei];
201 
202  // Get the list of labels
203  const labelList& curPoints = curFace;
204 
205  // For all points
206  forAll(curPoints, pointi)
207  {
208  // dGget the list of cells sharing this point
209  const labelList& curNeighbours =
210  PointCells[curPoints[pointi]];
211 
212  // For all neighbours
213  forAll(curNeighbours, neiI)
214  {
215  label curNei = curNeighbours[neiI];
216 
217  // Reject neighbours with the lower label
218  if (curNei > celli)
219  {
220  // Get the list of search faces
221  const faceList& searchFaces = cellsFaceShapes[curNei];
222 
223  forAll(searchFaces, neiFacei)
224  {
225  if (searchFaces[neiFacei] == curFace)
226  {
227  // Match!!
228  found = true;
229 
230  // Record the neighbour cell and face
231  neiCells[facei] = curNei;
232  faceOfNeiCell[facei] = neiFacei;
233  nNeighbours++;
234 
235  break;
236  }
237  }
238  if (found) break;
239  }
240  if (found) break;
241  }
242  if (found) break;
243  } // End of current points
244  } // End of current faces
245 
246  // Add the faces in the increasing order of neighbours
247  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
248  {
249  // Find the lowest neighbour which is still valid
250  label nextNei = -1;
251  label minNei = cells.size();
252 
253  forAll(neiCells, ncI)
254  {
255  if (neiCells[ncI] > -1 && neiCells[ncI] < minNei)
256  {
257  nextNei = ncI;
258  minNei = neiCells[ncI];
259  }
260  }
261 
262  if (nextNei > -1)
263  {
264  // Add the face to the list of faces
265  faces_[nFaces] = curFaces[nextNei];
266 
267  // Set cell-face and cell-neighbour-face to current face label
268  cells[celli][nextNei] = nFaces;
269  cells[neiCells[nextNei]][faceOfNeiCell[nextNei]] = nFaces;
270 
271  // Stop the neighbour from being used again
272  neiCells[nextNei] = -1;
273 
274  // Increment number of faces counter
275  nFaces++;
276  }
277  else
278  {
280  << "Error in internal face insertion"
281  << abort(FatalError);
282  }
283  }
284  }
285 
286  // Do boundary faces
287 
288  patchSizes.setSize(boundaryFaces.size(), -1);
289  patchStarts.setSize(boundaryFaces.size(), -1);
290 
291  forAll(boundaryFaces, patchi)
292  {
293  const faceList& patchFaces = boundaryFaces[patchi];
294 
295  labelList curPatchFaceCells =
296  facePatchFaceCells
297  (
298  patchFaces,
299  PointCells,
300  cellsFaceShapes,
301  patchi
302  );
303 
304  // Grab the start label
305  label curPatchStart = nFaces;
306 
307  forAll(patchFaces, facei)
308  {
309  const face& curFace = patchFaces[facei];
310 
311  const label cellInside = curPatchFaceCells[facei];
312 
313  // Get faces of the cell inside
314  const faceList& facesOfCellInside = cellsFaceShapes[cellInside];
315 
316  bool found = false;
317 
318  forAll(facesOfCellInside, cellFacei)
319  {
320  if (face::sameVertices(facesOfCellInside[cellFacei], curFace))
321  {
322  if (cells[cellInside][cellFacei] >= 0)
323  {
325  << "Trying to specify a boundary face " << curFace
326  << " on the face on cell " << cellInside
327  << " which is either an internal face or already "
328  << "belongs to some other patch. This is face "
329  << facei << " of patch "
330  << patchi << " named "
331  << boundaryPatchNames[patchi] << "."
332  << abort(FatalError);
333  }
334 
335  found = true;
336 
337  // Set the patch face to corresponding cell-face
338  faces_[nFaces] = facesOfCellInside[cellFacei];
339 
340  cells[cellInside][cellFacei] = nFaces;
341 
342  break;
343  }
344  }
345 
346  if (!found)
347  {
349  << "face " << facei << " of patch " << patchi
350  << " does not seem to belong to cell " << cellInside
351  << " which, according to the addressing, "
352  << "should be next to it."
353  << abort(FatalError);
354  }
355 
356  // Increment the counter of faces
357  nFaces++;
358  }
359 
360  patchSizes[patchi] = nFaces - curPatchStart;
361  patchStarts[patchi] = curPatchStart;
362  }
363 
364  // Grab "non-existing" faces and put them into a default patch
365 
366  defaultPatchStart = nFaces;
367 
368  forAll(cells, celli)
369  {
370  labelList& curCellFaces = cells[celli];
371 
372  forAll(curCellFaces, facei)
373  {
374  if (curCellFaces[facei] == -1) // "non-existent" face
375  {
376  curCellFaces[facei] = nFaces;
377  faces_[nFaces] = cellsFaceShapes[celli][facei];
378 
379  nFaces++;
380  }
381  }
382  }
383 
384  // Reset the size of the face list
385  faces_.setSize(nFaces);
386 }
387 
388 
390 (
391  const IOobject& io,
392  pointField&& points,
393  const cellShapeList& cellsAsShapes,
394  const faceListList& boundaryFaces,
395  const wordList& boundaryPatchNames,
396  const wordList& boundaryPatchTypes,
397  const word& defaultBoundaryPatchName,
398  const word& defaultBoundaryPatchType,
399  const wordList& boundaryPatchPhysicalTypes,
400  const bool syncPar
401 )
402 :
403  objectRegistry(io),
404  primitiveMesh(),
405  points_
406  (
407  IOobject
408  (
409  "points",
410  instance(),
411  meshSubDir,
412  *this,
413  IOobject::NO_READ,
414  IOobject::AUTO_WRITE
415  ),
416  move(points)
417  ),
418  faces_
419  (
420  IOobject
421  (
422  "faces",
423  instance(),
424  meshSubDir,
425  *this,
426  IOobject::NO_READ,
427  IOobject::AUTO_WRITE
428  ),
429  label(0)
430  ),
431  owner_
432  (
433  IOobject
434  (
435  "owner",
436  instance(),
437  meshSubDir,
438  *this,
439  IOobject::NO_READ,
440  IOobject::AUTO_WRITE
441  ),
442  label(0)
443  ),
444  neighbour_
445  (
446  IOobject
447  (
448  "neighbour",
449  instance(),
450  meshSubDir,
451  *this,
452  IOobject::NO_READ,
453  IOobject::AUTO_WRITE
454  ),
455  label(0)
456  ),
457  clearedPrimitives_(false),
458  boundary_
459  (
460  IOobject
461  (
462  "boundary",
463  instance(),
464  meshSubDir,
465  *this,
466  IOobject::NO_READ,
467  IOobject::AUTO_WRITE
468  ),
469  *this,
470  boundaryFaces.size() + 1 // Add room for a default patch
471  ),
472  bounds_(points_, syncPar),
473  comm_(UPstream::worldComm),
474  geometricD_(Zero),
475  solutionD_(Zero),
476  tetBasePtIsPtr_(nullptr),
477  cellTreePtr_(nullptr),
478  pointZones_
479  (
480  IOobject
481  (
482  "pointZones",
483  instance(),
484  meshSubDir,
485  *this,
486  IOobject::NO_READ,
487  IOobject::NO_WRITE
488  ),
489  *this
490  ),
491  faceZones_
492  (
493  IOobject
494  (
495  "faceZones",
496  instance(),
497  meshSubDir,
498  *this,
499  IOobject::NO_READ,
500  IOobject::NO_WRITE
501  ),
502  *this
503  ),
504  cellZones_
505  (
506  IOobject
507  (
508  "cellZones",
509  instance(),
510  meshSubDir,
511  *this,
512  IOobject::NO_READ,
513  IOobject::NO_WRITE
514  ),
515  *this
516  ),
517  globalMeshDataPtr_(nullptr),
518  curMotionTimeIndex_(-1),
519  oldPointsPtr_(nullptr),
520  oldCellCentresPtr_(nullptr),
521  storeOldCellCentres_(false),
522  moving_(false),
523  topoChanged_(false)
524 {
525  if (debug)
526  {
527  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
528  }
529 
530  // Remove all of the old mesh files if they exist
532 
533  // Calculate faces and cells
534  labelList patchSizes;
535  labelList patchStarts;
536  label defaultPatchStart;
537  label nFaces;
538  cellList cells;
539  setTopology
540  (
541  cellsAsShapes,
542  boundaryFaces,
543  boundaryPatchNames,
544  patchSizes,
545  patchStarts,
546  defaultPatchStart,
547  nFaces,
548  cells
549  );
550 
551  // Warning: Patches can only be added once the face list is
552  // completed, as they hold a subList of the face list
553  forAll(boundaryFaces, patchi)
554  {
555  // Add the patch to the list
556  boundary_.set
557  (
558  patchi,
560  (
561  boundaryPatchTypes[patchi],
562  boundaryPatchNames[patchi],
563  patchSizes[patchi],
564  patchStarts[patchi],
565  patchi,
566  boundary_
567  )
568  );
569 
570  if
571  (
572  boundaryPatchPhysicalTypes.size()
573  && boundaryPatchPhysicalTypes[patchi].size()
574  )
575  {
576  boundary_[patchi].physicalType() =
577  boundaryPatchPhysicalTypes[patchi];
578  }
579  }
580 
581  label nAllPatches = boundaryFaces.size();
582 
583 
584  label nDefaultFaces = nFaces - defaultPatchStart;
585  if (syncPar)
586  {
587  reduce(nDefaultFaces, sumOp<label>());
588  }
589 
590  if (nDefaultFaces > 0)
591  {
592  if (debug)
593  {
595  << "Found " << nDefaultFaces
596  << " undefined faces in mesh; adding to default patch." << endl;
597  }
598 
599  // Check if there already exists a defaultFaces patch as last patch
600  // and reuse it.
601  label patchi = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
602 
603  if (patchi != -1)
604  {
605  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
606  {
608  << "Default patch " << boundary_[patchi].name()
609  << " already has faces in it or is not"
610  << " last in list of patches." << exit(FatalError);
611  }
612 
614  << "Reusing existing patch " << patchi
615  << " for undefined faces." << endl;
616 
617  boundary_.set
618  (
619  patchi,
621  (
622  boundaryPatchTypes[patchi],
623  boundaryPatchNames[patchi],
624  nFaces - defaultPatchStart,
625  defaultPatchStart,
626  patchi,
627  boundary_
628  )
629  );
630  }
631  else
632  {
633  boundary_.set
634  (
635  nAllPatches,
637  (
638  defaultBoundaryPatchType,
639  defaultBoundaryPatchName,
640  nFaces - defaultPatchStart,
641  defaultPatchStart,
642  boundary_.size() - 1,
643  boundary_
644  )
645  );
646 
647  nAllPatches++;
648  }
649  }
650 
651  // Reset the size of the boundary
652  boundary_.setSize(nAllPatches);
653 
654  // Set the primitive mesh
655  initMesh(cells);
656 
657  if (syncPar)
658  {
659  // Calculate topology for the patches (processor-processor comms etc.)
660  boundary_.topoChange();
661 
662  // Calculate the geometry for the patches (transformation tensors etc.)
663  boundary_.calcGeometry();
664  }
665 }
666 
667 
669 (
670  const IOobject& io,
671  pointField&& points,
672  const cellShapeList& cellsAsShapes,
673  const faceListList& boundaryFaces,
674  const wordList& boundaryPatchNames,
675  const PtrList<dictionary>& boundaryDicts,
676  const word& defaultBoundaryPatchName,
677  const word& defaultBoundaryPatchType,
678  const bool syncPar
679 )
680 :
681  objectRegistry(io),
682  primitiveMesh(),
683  points_
684  (
685  IOobject
686  (
687  "points",
688  instance(),
689  meshSubDir,
690  *this,
691  IOobject::NO_READ,
692  IOobject::AUTO_WRITE
693  ),
694  move(points)
695  ),
696  faces_
697  (
698  IOobject
699  (
700  "faces",
701  instance(),
702  meshSubDir,
703  *this,
704  IOobject::NO_READ,
705  IOobject::AUTO_WRITE
706  ),
707  label(0)
708  ),
709  owner_
710  (
711  IOobject
712  (
713  "owner",
714  instance(),
715  meshSubDir,
716  *this,
717  IOobject::NO_READ,
718  IOobject::AUTO_WRITE
719  ),
720  label(0)
721  ),
722  neighbour_
723  (
724  IOobject
725  (
726  "neighbour",
727  instance(),
728  meshSubDir,
729  *this,
730  IOobject::NO_READ,
731  IOobject::AUTO_WRITE
732  ),
733  label(0)
734  ),
735  clearedPrimitives_(false),
736  boundary_
737  (
738  IOobject
739  (
740  "boundary",
741  instance(),
742  meshSubDir,
743  *this,
744  IOobject::NO_READ,
745  IOobject::AUTO_WRITE
746  ),
747  *this,
748  boundaryFaces.size() + 1 // Add room for a default patch
749  ),
750  bounds_(points_, syncPar),
751  comm_(UPstream::worldComm),
752  geometricD_(Zero),
753  solutionD_(Zero),
754  tetBasePtIsPtr_(nullptr),
755  cellTreePtr_(nullptr),
756  pointZones_
757  (
758  IOobject
759  (
760  "pointZones",
761  instance(),
762  meshSubDir,
763  *this,
764  IOobject::NO_READ,
765  IOobject::NO_WRITE
766  ),
767  *this
768  ),
769  faceZones_
770  (
771  IOobject
772  (
773  "faceZones",
774  instance(),
775  meshSubDir,
776  *this,
777  IOobject::NO_READ,
778  IOobject::NO_WRITE
779  ),
780  *this
781  ),
782  cellZones_
783  (
784  IOobject
785  (
786  "cellZones",
787  instance(),
788  meshSubDir,
789  *this,
790  IOobject::NO_READ,
791  IOobject::NO_WRITE
792  ),
793  *this
794  ),
795  globalMeshDataPtr_(nullptr),
796  curMotionTimeIndex_(-1),
797  oldPointsPtr_(nullptr),
798  oldCellCentresPtr_(nullptr),
799  storeOldCellCentres_(false),
800  moving_(false)
801 {
802  if (debug)
803  {
804  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
805  }
806 
807  // Remove all of the old mesh files if they exist
809 
810  // Calculate faces and cells
811  labelList patchSizes;
812  labelList patchStarts;
813  label defaultPatchStart;
814  label nFaces;
815  cellList cells;
816  setTopology
817  (
818  cellsAsShapes,
819  boundaryFaces,
820  boundaryPatchNames,
821  patchSizes,
822  patchStarts,
823  defaultPatchStart,
824  nFaces,
825  cells
826  );
827 
828  // Warning: Patches can only be added once the face list is
829  // completed, as they hold a subList of the face list
830  forAll(boundaryDicts, patchi)
831  {
832  dictionary patchDict(boundaryDicts[patchi]);
833 
834  patchDict.set("nFaces", patchSizes[patchi]);
835  patchDict.set("startFace", patchStarts[patchi]);
836 
837  // Add the patch to the list
838  boundary_.set
839  (
840  patchi,
842  (
843  boundaryPatchNames[patchi],
844  patchDict,
845  patchi,
846  boundary_
847  )
848  );
849  }
850 
851  label nAllPatches = boundaryFaces.size();
852 
853  label nDefaultFaces = nFaces - defaultPatchStart;
854  if (syncPar)
855  {
856  reduce(nDefaultFaces, sumOp<label>());
857  }
858 
859  if (nDefaultFaces > 0)
860  {
861  if (debug)
862  {
864  << "Found " << nDefaultFaces
865  << " undefined faces in mesh; adding to default patch." << endl;
866  }
867 
868  // Check if there already exists a defaultFaces patch as last patch
869  // and reuse it.
870  label patchi = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
871 
872  if (patchi != -1)
873  {
874  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
875  {
877  << "Default patch " << boundary_[patchi].name()
878  << " already has faces in it or is not"
879  << " last in list of patches." << exit(FatalError);
880  }
881 
883  << "Reusing existing patch " << patchi
884  << " for undefined faces." << endl;
885 
886  boundary_.set
887  (
888  patchi,
890  (
891  boundary_[patchi].type(),
892  boundary_[patchi].name(),
893  nFaces - defaultPatchStart,
894  defaultPatchStart,
895  patchi,
896  boundary_
897  )
898  );
899  }
900  else
901  {
902  boundary_.set
903  (
904  nAllPatches,
906  (
907  defaultBoundaryPatchType,
908  defaultBoundaryPatchName,
909  nFaces - defaultPatchStart,
910  defaultPatchStart,
911  boundary_.size() - 1,
912  boundary_
913  )
914  );
915 
916  nAllPatches++;
917  }
918  }
919 
920  // Reset the size of the boundary
921  boundary_.setSize(nAllPatches);
922 
923  // Set the primitive mesh
924  initMesh(cells);
925 
926  if (syncPar)
927  {
928  // Calculate topology for the patches (processor-processor comms etc.)
929  boundary_.topoChange();
930 
931  // Calculate the geometry for the patches (transformation tensors etc.)
932  boundary_.calcGeometry();
933  }
934 }
935 
936 
937 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Inter-processor communications stream.
Definition: UPstream.H:59
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1152
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:150
Registry of regIOobjects.
void topoChange()
Correct polyBoundaryMesh after topology update.
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:162
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1569
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
label nFaces() const
const cellList & cells() const
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
List< faceList > faceListList
Definition: faceListFwd.H:43
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
List< face > faceList
Definition: faceListFwd.H:41
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488