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-2022 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,
415  ),
416  move(points)
417  ),
418  faces_
419  (
420  IOobject
421  (
422  "faces",
423  instance(),
424  meshSubDir,
425  *this,
428  ),
429  0
430  ),
431  owner_
432  (
433  IOobject
434  (
435  "owner",
436  instance(),
437  meshSubDir,
438  *this,
441  ),
442  0
443  ),
444  neighbour_
445  (
446  IOobject
447  (
448  "neighbour",
449  instance(),
450  meshSubDir,
451  *this,
454  ),
455  0
456  ),
457  clearedPrimitives_(false),
458  boundary_
459  (
460  IOobject
461  (
462  "boundary",
463  instance(),
464  meshSubDir,
465  *this,
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,
488  ),
489  *this,
490  0
491  ),
492  faceZones_
493  (
494  IOobject
495  (
496  "faceZones",
497  instance(),
498  meshSubDir,
499  *this,
502  ),
503  *this,
504  0
505  ),
506  cellZones_
507  (
508  IOobject
509  (
510  "cellZones",
511  instance(),
512  meshSubDir,
513  *this,
516  ),
517  *this,
518  0
519  ),
520  globalMeshDataPtr_(nullptr),
521  curMotionTimeIndex_(-1),
522  oldPointsPtr_(nullptr),
523  oldCellCentresPtr_(nullptr),
524  storeOldCellCentres_(false),
525  moving_(false),
526  topoChanged_(false)
527 {
528  if (debug)
529  {
530  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
531  }
532 
533  // Remove all of the old mesh files if they exist
535 
536  // Calculate faces and cells
537  labelList patchSizes;
538  labelList patchStarts;
539  label defaultPatchStart;
540  label nFaces;
541  cellList cells;
542  setTopology
543  (
544  cellsAsShapes,
545  boundaryFaces,
546  boundaryPatchNames,
547  patchSizes,
548  patchStarts,
549  defaultPatchStart,
550  nFaces,
551  cells
552  );
553 
554  // Warning: Patches can only be added once the face list is
555  // completed, as they hold a subList of the face list
556  forAll(boundaryFaces, patchi)
557  {
558  // Add the patch to the list
559  boundary_.set
560  (
561  patchi,
563  (
564  boundaryPatchTypes[patchi],
565  boundaryPatchNames[patchi],
566  patchSizes[patchi],
567  patchStarts[patchi],
568  patchi,
569  boundary_
570  )
571  );
572 
573  if
574  (
575  boundaryPatchPhysicalTypes.size()
576  && boundaryPatchPhysicalTypes[patchi].size()
577  )
578  {
579  boundary_[patchi].physicalType() =
580  boundaryPatchPhysicalTypes[patchi];
581  }
582  }
583 
584  label nAllPatches = boundaryFaces.size();
585 
586 
587  label nDefaultFaces = nFaces - defaultPatchStart;
588  if (syncPar)
589  {
590  reduce(nDefaultFaces, sumOp<label>());
591  }
592 
593  if (nDefaultFaces > 0)
594  {
595  if (debug)
596  {
598  << "Found " << nDefaultFaces
599  << " undefined faces in mesh; adding to default patch." << endl;
600  }
601 
602  // Check if there already exists a defaultFaces patch as last patch
603  // and reuse it.
604  label patchi = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
605 
606  if (patchi != -1)
607  {
608  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
609  {
611  << "Default patch " << boundary_[patchi].name()
612  << " already has faces in it or is not"
613  << " last in list of patches." << exit(FatalError);
614  }
615 
617  << "Reusing existing patch " << patchi
618  << " for undefined faces." << endl;
619 
620  boundary_.set
621  (
622  patchi,
624  (
625  boundaryPatchTypes[patchi],
626  boundaryPatchNames[patchi],
627  nFaces - defaultPatchStart,
628  defaultPatchStart,
629  patchi,
630  boundary_
631  )
632  );
633  }
634  else
635  {
636  boundary_.set
637  (
638  nAllPatches,
640  (
641  defaultBoundaryPatchType,
642  defaultBoundaryPatchName,
643  nFaces - defaultPatchStart,
644  defaultPatchStart,
645  boundary_.size() - 1,
646  boundary_
647  )
648  );
649 
650  nAllPatches++;
651  }
652  }
653 
654  // Reset the size of the boundary
655  boundary_.setSize(nAllPatches);
656 
657  // Set the primitive mesh
658  initMesh(cells);
659 
660  if (syncPar)
661  {
662  // Calculate topology for the patches (processor-processor comms etc.)
663  boundary_.topoChange();
664 
665  // Calculate the geometry for the patches (transformation tensors etc.)
666  boundary_.calcGeometry();
667  }
668 
669  if (debug)
670  {
671  if (checkMesh())
672  {
673  Info<< "Mesh OK" << endl;
674  }
675  }
676 }
677 
678 
680 (
681  const IOobject& io,
682  pointField&& points,
683  const cellShapeList& cellsAsShapes,
684  const faceListList& boundaryFaces,
685  const wordList& boundaryPatchNames,
686  const PtrList<dictionary>& boundaryDicts,
687  const word& defaultBoundaryPatchName,
688  const word& defaultBoundaryPatchType,
689  const bool syncPar
690 )
691 :
692  objectRegistry(io),
693  primitiveMesh(),
694  points_
695  (
696  IOobject
697  (
698  "points",
699  instance(),
700  meshSubDir,
701  *this,
704  ),
705  move(points)
706  ),
707  faces_
708  (
709  IOobject
710  (
711  "faces",
712  instance(),
713  meshSubDir,
714  *this,
717  ),
718  0
719  ),
720  owner_
721  (
722  IOobject
723  (
724  "owner",
725  instance(),
726  meshSubDir,
727  *this,
730  ),
731  0
732  ),
733  neighbour_
734  (
735  IOobject
736  (
737  "neighbour",
738  instance(),
739  meshSubDir,
740  *this,
743  ),
744  0
745  ),
746  clearedPrimitives_(false),
747  boundary_
748  (
749  IOobject
750  (
751  "boundary",
752  instance(),
753  meshSubDir,
754  *this,
757  ),
758  *this,
759  boundaryFaces.size() + 1 // Add room for a default patch
760  ),
761  bounds_(points_, syncPar),
762  comm_(UPstream::worldComm),
763  geometricD_(Zero),
764  solutionD_(Zero),
765  tetBasePtIsPtr_(nullptr),
766  cellTreePtr_(nullptr),
767  pointZones_
768  (
769  IOobject
770  (
771  "pointZones",
772  instance(),
773  meshSubDir,
774  *this,
777  ),
778  *this,
779  0
780  ),
781  faceZones_
782  (
783  IOobject
784  (
785  "faceZones",
786  instance(),
787  meshSubDir,
788  *this,
791  ),
792  *this,
793  0
794  ),
795  cellZones_
796  (
797  IOobject
798  (
799  "cellZones",
800  instance(),
801  meshSubDir,
802  *this,
805  ),
806  *this,
807  0
808  ),
809  globalMeshDataPtr_(nullptr),
810  curMotionTimeIndex_(-1),
811  oldPointsPtr_(nullptr),
812  oldCellCentresPtr_(nullptr),
813  storeOldCellCentres_(false),
814  moving_(false),
815  topoChanged_(false)
816 {
817  if (debug)
818  {
819  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
820  }
821 
822  // Remove all of the old mesh files if they exist
824 
825  // Calculate faces and cells
826  labelList patchSizes;
827  labelList patchStarts;
828  label defaultPatchStart;
829  label nFaces;
830  cellList cells;
831  setTopology
832  (
833  cellsAsShapes,
834  boundaryFaces,
835  boundaryPatchNames,
836  patchSizes,
837  patchStarts,
838  defaultPatchStart,
839  nFaces,
840  cells
841  );
842 
843  // Warning: Patches can only be added once the face list is
844  // completed, as they hold a subList of the face list
845  forAll(boundaryDicts, patchi)
846  {
847  dictionary patchDict(boundaryDicts[patchi]);
848 
849  patchDict.set("nFaces", patchSizes[patchi]);
850  patchDict.set("startFace", patchStarts[patchi]);
851 
852  // Add the patch to the list
853  boundary_.set
854  (
855  patchi,
857  (
858  boundaryPatchNames[patchi],
859  patchDict,
860  patchi,
861  boundary_
862  )
863  );
864  }
865 
866  label nAllPatches = boundaryFaces.size();
867 
868  label nDefaultFaces = nFaces - defaultPatchStart;
869  if (syncPar)
870  {
871  reduce(nDefaultFaces, sumOp<label>());
872  }
873 
874  if (nDefaultFaces > 0)
875  {
876  if (debug)
877  {
879  << "Found " << nDefaultFaces
880  << " undefined faces in mesh; adding to default patch." << endl;
881  }
882 
883  // Check if there already exists a defaultFaces patch as last patch
884  // and reuse it.
885  label patchi = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
886 
887  if (patchi != -1)
888  {
889  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
890  {
892  << "Default patch " << boundary_[patchi].name()
893  << " already has faces in it or is not"
894  << " last in list of patches." << exit(FatalError);
895  }
896 
898  << "Reusing existing patch " << patchi
899  << " for undefined faces." << endl;
900 
901  boundary_.set
902  (
903  patchi,
905  (
906  boundary_[patchi].type(),
907  boundary_[patchi].name(),
908  nFaces - defaultPatchStart,
909  defaultPatchStart,
910  patchi,
911  boundary_
912  )
913  );
914  }
915  else
916  {
917  boundary_.set
918  (
919  nAllPatches,
921  (
922  defaultBoundaryPatchType,
923  defaultBoundaryPatchName,
924  nFaces - defaultPatchStart,
925  defaultPatchStart,
926  boundary_.size() - 1,
927  boundary_
928  )
929  );
930 
931  nAllPatches++;
932  }
933  }
934 
935  // Reset the size of the boundary
936  boundary_.setSize(nAllPatches);
937 
938  // Set the primitive mesh
939  initMesh(cells);
940 
941  if (syncPar)
942  {
943  // Calculate topology for the patches (processor-processor comms etc.)
944  boundary_.topoChange();
945 
946  // Calculate the geometry for the patches (transformation tensors etc.)
947  boundary_.calcGeometry();
948  }
949 
950  if (debug)
951  {
952  if (checkMesh())
953  {
954  Info << "Mesh OK" << endl;
955  }
956  }
957 }
958 
959 
960 // ************************************************************************* //
bool moving_
Member data pending transfer to fvMesh.
Definition: polyMesh.H:307
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
List< faceList > faceListList
Definition: faceListFwd.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
label nFaces() const
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
primitiveMesh()
Construct null.
Definition: primitiveMesh.C:39
const cellList & cells() const
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:278
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
objectRegistry(const Time &db, const label nIoObjects=128)
Construct the time objectRegistry given an initial estimate.
bool topoChanged_
Is the mesh topology changing.
Definition: polyMesh.H:310
void topoChange()
Correct polyBoundaryMesh after topology update.
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:43
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1489
bool found(const word &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
virtual bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static bool sameVertices(const face &, const face &)
Return true if the faces have the same vertices.
Definition: face.C:150
List< word > wordList
A List of words.
Definition: fileName.H:54
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1291
messageStream Info
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
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
List< cell > cellList
list of cells
Definition: cellList.H:42
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163