polyMeshFromShapeMesh.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 \*---------------------------------------------------------------------------*/
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 
37 (
38  const cellShapeList& c
39 ) const
40 {
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];
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 
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 //- Set faces_, calculate cells and patchStarts.
132 (
133  const cellShapeList& cellsAsShapes,
134  const faceListList& boundaryFaces,
135  const wordList& boundaryPatchNames,
136  labelList& patchSizes,
137  labelList& patchStarts,
138  label& defaultPatchStart,
139  label& nFaces,
140  cellList& cells
141 )
142 {
143  // Calculate the faces of all cells
144  // Initialise maximum possible numer of mesh faces to 0
145  label maxFaces = 0;
146 
147  // Set up a list of face shapes for each cell
148  faceListList cellsFaceShapes(cellsAsShapes.size());
149  cells.setSize(cellsAsShapes.size());
150 
151  forAll(cellsFaceShapes, cellI)
152  {
153  cellsFaceShapes[cellI] = cellsAsShapes[cellI].faces();
154 
155  cells[cellI].setSize(cellsFaceShapes[cellI].size());
156 
157  // Initialise cells to -1 to flag undefined faces
158  static_cast<labelList&>(cells[cellI]) = -1;
159 
160  // Count maximum possible numer of mesh faces
161  maxFaces += cellsFaceShapes[cellI].size();
162  }
163 
164  // Set size of faces array to maximum possible number of mesh faces
165  faces_.setSize(maxFaces);
166 
167  // Initialise number of faces to 0
168  nFaces = 0;
169 
170  // Set reference to point-cell addressing
171  labelListList PointCells = cellShapePointCells(cellsAsShapes);
172 
173  bool found = false;
174 
175  forAll(cells, cellI)
176  {
177  // Note:
178  // Insertion cannot be done in one go as the faces need to be
179  // added into the list in the increasing order of neighbour
180  // cells. Therefore, all neighbours will be detected first
181  // and then added in the correct order.
182 
183  const faceList& curFaces = cellsFaceShapes[cellI];
184 
185  // Record the neighbour cell
186  labelList neiCells(curFaces.size(), -1);
187 
188  // Record the face of neighbour cell
189  labelList faceOfNeiCell(curFaces.size(), -1);
190 
191  label nNeighbours = 0;
192 
193  // For all faces ...
194  forAll(curFaces, faceI)
195  {
196  // Skip faces that have already been matched
197  if (cells[cellI][faceI] >= 0) continue;
198 
199  found = false;
200 
201  const face& curFace = curFaces[faceI];
202 
203  // Get the list of labels
204  const labelList& curPoints = curFace;
205 
206  // For all points
207  forAll(curPoints, pointI)
208  {
209  // dGget the list of cells sharing this point
210  const labelList& curNeighbours =
211  PointCells[curPoints[pointI]];
212 
213  // For all neighbours
214  forAll(curNeighbours, neiI)
215  {
216  label curNei = curNeighbours[neiI];
217 
218  // Reject neighbours with the lower label
219  if (curNei > cellI)
220  {
221  // Get the list of search faces
222  const faceList& searchFaces = cellsFaceShapes[curNei];
223 
224  forAll(searchFaces, neiFaceI)
225  {
226  if (searchFaces[neiFaceI] == curFace)
227  {
228  // Match!!
229  found = true;
230 
231  // Record the neighbour cell and face
232  neiCells[faceI] = curNei;
233  faceOfNeiCell[faceI] = neiFaceI;
234  nNeighbours++;
235 
236  break;
237  }
238  }
239  if (found) break;
240  }
241  if (found) break;
242  }
243  if (found) break;
244  } // End of current points
245  } // End of current faces
246 
247  // Add the faces in the increasing order of neighbours
248  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
249  {
250  // Find the lowest neighbour which is still valid
251  label nextNei = -1;
252  label minNei = cells.size();
253 
254  forAll(neiCells, ncI)
255  {
256  if (neiCells[ncI] > -1 && neiCells[ncI] < minNei)
257  {
258  nextNei = ncI;
259  minNei = neiCells[ncI];
260  }
261  }
262 
263  if (nextNei > -1)
264  {
265  // Add the face to the list of faces
266  faces_[nFaces] = curFaces[nextNei];
267 
268  // Set cell-face and cell-neighbour-face to current face label
269  cells[cellI][nextNei] = nFaces;
270  cells[neiCells[nextNei]][faceOfNeiCell[nextNei]] = nFaces;
271 
272  // Stop the neighbour from being used again
273  neiCells[nextNei] = -1;
274 
275  // Increment number of faces counter
276  nFaces++;
277  }
278  else
279  {
281  << "Error in internal face insertion"
282  << abort(FatalError);
283  }
284  }
285  }
286 
287  // Do boundary faces
288 
289  patchSizes.setSize(boundaryFaces.size(), -1);
290  patchStarts.setSize(boundaryFaces.size(), -1);
291 
292  forAll(boundaryFaces, patchI)
293  {
294  const faceList& patchFaces = boundaryFaces[patchI];
295 
296  labelList curPatchFaceCells =
297  facePatchFaceCells
298  (
299  patchFaces,
300  PointCells,
301  cellsFaceShapes,
302  patchI
303  );
304 
305  // Grab the start label
306  label curPatchStart = nFaces;
307 
308  forAll(patchFaces, faceI)
309  {
310  const face& curFace = patchFaces[faceI];
311 
312  const label cellInside = curPatchFaceCells[faceI];
313 
314  // Get faces of the cell inside
315  const faceList& facesOfCellInside = cellsFaceShapes[cellInside];
316 
317  bool found = false;
318 
319  forAll(facesOfCellInside, cellFaceI)
320  {
321  if (face::sameVertices(facesOfCellInside[cellFaceI], curFace))
322  {
323  if (cells[cellInside][cellFaceI] >= 0)
324  {
326  << "Trying to specify a boundary face " << curFace
327  << " on the face on cell " << cellInside
328  << " which is either an internal face or already "
329  << "belongs to some other patch. This is face "
330  << faceI << " of patch "
331  << patchI << " named "
332  << boundaryPatchNames[patchI] << "."
333  << abort(FatalError);
334  }
335 
336  found = true;
337 
338  // Set the patch face to corresponding cell-face
339  faces_[nFaces] = facesOfCellInside[cellFaceI];
340 
341  cells[cellInside][cellFaceI] = nFaces;
342 
343  break;
344  }
345  }
346 
347  if (!found)
348  {
350  << "face " << faceI << " of patch " << patchI
351  << " does not seem to belong to cell " << cellInside
352  << " which, according to the addressing, "
353  << "should be next to it."
354  << abort(FatalError);
355  }
356 
357  // Increment the counter of faces
358  nFaces++;
359  }
360 
361  patchSizes[patchI] = nFaces - curPatchStart;
362  patchStarts[patchI] = curPatchStart;
363  }
364 
365  // Grab "non-existing" faces and put them into a default patch
366 
367  defaultPatchStart = nFaces;
368 
369  forAll(cells, cellI)
370  {
371  labelList& curCellFaces = cells[cellI];
372 
373  forAll(curCellFaces, faceI)
374  {
375  if (curCellFaces[faceI] == -1) // "non-existent" face
376  {
377  curCellFaces[faceI] = nFaces;
378  faces_[nFaces] = cellsFaceShapes[cellI][faceI];
379 
380  nFaces++;
381  }
382  }
383  }
384 
385  // Reset the size of the face list
386  faces_.setSize(nFaces);
387 
388  return ;
389 }
390 
391 
393 (
394  const IOobject& io,
395  const Xfer<pointField>& points,
396  const cellShapeList& cellsAsShapes,
397  const faceListList& boundaryFaces,
398  const wordList& boundaryPatchNames,
399  const wordList& boundaryPatchTypes,
400  const word& defaultBoundaryPatchName,
401  const word& defaultBoundaryPatchType,
402  const wordList& boundaryPatchPhysicalTypes,
403  const bool syncPar
404 )
405 :
406  objectRegistry(io),
407  primitiveMesh(),
408  points_
409  (
410  IOobject
411  (
412  "points",
413  instance(),
414  meshSubDir,
415  *this,
416  IOobject::NO_READ,
417  IOobject::AUTO_WRITE
418  ),
419  points
420  ),
421  faces_
422  (
423  IOobject
424  (
425  "faces",
426  instance(),
427  meshSubDir,
428  *this,
429  IOobject::NO_READ,
430  IOobject::AUTO_WRITE
431  ),
432  0
433  ),
434  owner_
435  (
436  IOobject
437  (
438  "owner",
439  instance(),
440  meshSubDir,
441  *this,
442  IOobject::NO_READ,
443  IOobject::AUTO_WRITE
444  ),
445  0
446  ),
447  neighbour_
448  (
449  IOobject
450  (
451  "neighbour",
452  instance(),
453  meshSubDir,
454  *this,
455  IOobject::NO_READ,
456  IOobject::AUTO_WRITE
457  ),
458  0
459  ),
460  clearedPrimitives_(false),
461  boundary_
462  (
463  IOobject
464  (
465  "boundary",
466  instance(),
467  meshSubDir,
468  *this,
469  IOobject::NO_READ,
470  IOobject::AUTO_WRITE
471  ),
472  *this,
473  boundaryFaces.size() + 1 // Add room for a default patch
474  ),
475  bounds_(points_, syncPar),
476  comm_(UPstream::worldComm),
477  geometricD_(Vector<label>::zero),
478  solutionD_(Vector<label>::zero),
479  tetBasePtIsPtr_(NULL),
480  cellTreePtr_(NULL),
481  pointZones_
482  (
483  IOobject
484  (
485  "pointZones",
486  instance(),
487  meshSubDir,
488  *this,
489  IOobject::NO_READ,
490  IOobject::NO_WRITE
491  ),
492  *this,
493  0
494  ),
495  faceZones_
496  (
497  IOobject
498  (
499  "faceZones",
500  instance(),
501  meshSubDir,
502  *this,
503  IOobject::NO_READ,
504  IOobject::NO_WRITE
505  ),
506  *this,
507  0
508  ),
509  cellZones_
510  (
511  IOobject
512  (
513  "cellZones",
514  instance(),
515  meshSubDir,
516  *this,
517  IOobject::NO_READ,
518  IOobject::NO_WRITE
519  ),
520  *this,
521  0
522  ),
523  globalMeshDataPtr_(NULL),
524  moving_(false),
525  topoChanging_(false),
526  curMotionTimeIndex_(time().timeIndex()),
527  oldPointsPtr_(NULL)
528 {
529  if (debug)
530  {
531  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
532  }
533 
534  // Remove all of the old mesh files if they exist
535  removeFiles(instance());
536 
537  // Calculate faces and cells
538  labelList patchSizes;
539  labelList patchStarts;
540  label defaultPatchStart;
541  label nFaces;
542  cellList cells;
543  setTopology
544  (
545  cellsAsShapes,
546  boundaryFaces,
547  boundaryPatchNames,
548  patchSizes,
549  patchStarts,
550  defaultPatchStart,
551  nFaces,
552  cells
553  );
554 
555  // Warning: Patches can only be added once the face list is
556  // completed, as they hold a subList of the face list
557  forAll(boundaryFaces, patchI)
558  {
559  // Add the patch to the list
560  boundary_.set
561  (
562  patchI,
564  (
565  boundaryPatchTypes[patchI],
566  boundaryPatchNames[patchI],
567  patchSizes[patchI],
568  patchStarts[patchI],
569  patchI,
570  boundary_
571  )
572  );
573 
574  if
575  (
576  boundaryPatchPhysicalTypes.size()
577  && boundaryPatchPhysicalTypes[patchI].size()
578  )
579  {
580  boundary_[patchI].physicalType() =
581  boundaryPatchPhysicalTypes[patchI];
582  }
583  }
584 
585  label nAllPatches = boundaryFaces.size();
586 
587 
588  label nDefaultFaces = nFaces - defaultPatchStart;
589  if (syncPar)
590  {
591  reduce(nDefaultFaces, sumOp<label>());
592  }
593 
594  if (nDefaultFaces > 0)
595  {
597  << "Found " << nDefaultFaces
598  << " undefined faces in mesh; adding to default patch." << endl;
599 
600  // Check if there already exists a defaultFaces patch as last patch
601  // and reuse it.
602  label patchI = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
603 
604  if (patchI != -1)
605  {
606  if (patchI != boundaryFaces.size()-1 || boundary_[patchI].size())
607  {
609  << "Default patch " << boundary_[patchI].name()
610  << " already has faces in it or is not"
611  << " last in list of patches." << exit(FatalError);
612  }
613 
615  << "Reusing existing patch " << patchI
616  << " for undefined faces." << endl;
617 
618  boundary_.set
619  (
620  patchI,
622  (
623  boundaryPatchTypes[patchI],
624  boundaryPatchNames[patchI],
625  nFaces - defaultPatchStart,
626  defaultPatchStart,
627  patchI,
628  boundary_
629  )
630  );
631  }
632  else
633  {
634  boundary_.set
635  (
636  nAllPatches,
638  (
639  defaultBoundaryPatchType,
640  defaultBoundaryPatchName,
641  nFaces - defaultPatchStart,
642  defaultPatchStart,
643  boundary_.size() - 1,
644  boundary_
645  )
646  );
647 
648  nAllPatches++;
649  }
650  }
651 
652  // Reset the size of the boundary
653  boundary_.setSize(nAllPatches);
654 
655  // Set the primitive mesh
656  initMesh(cells);
657 
658  if (syncPar)
659  {
660  // Calculate topology for the patches (processor-processor comms etc.)
661  boundary_.updateMesh();
662 
663  // Calculate the geometry for the patches (transformation tensors etc.)
664  boundary_.calcGeometry();
665  }
666 
667  if (debug)
668  {
669  if (checkMesh())
670  {
671  Info<< "Mesh OK" << endl;
672  }
673  }
674 }
675 
676 
678 (
679  const IOobject& io,
680  const Xfer<pointField>& points,
681  const cellShapeList& cellsAsShapes,
682  const faceListList& boundaryFaces,
683  const wordList& boundaryPatchNames,
684  const PtrList<dictionary>& boundaryDicts,
685  const word& defaultBoundaryPatchName,
686  const word& defaultBoundaryPatchType,
687  const bool syncPar
688 )
689 :
690  objectRegistry(io),
691  primitiveMesh(),
692  points_
693  (
694  IOobject
695  (
696  "points",
697  instance(),
698  meshSubDir,
699  *this,
700  IOobject::NO_READ,
701  IOobject::AUTO_WRITE
702  ),
703  points
704  ),
705  faces_
706  (
707  IOobject
708  (
709  "faces",
710  instance(),
711  meshSubDir,
712  *this,
713  IOobject::NO_READ,
714  IOobject::AUTO_WRITE
715  ),
716  0
717  ),
718  owner_
719  (
720  IOobject
721  (
722  "owner",
723  instance(),
724  meshSubDir,
725  *this,
726  IOobject::NO_READ,
727  IOobject::AUTO_WRITE
728  ),
729  0
730  ),
731  neighbour_
732  (
733  IOobject
734  (
735  "neighbour",
736  instance(),
737  meshSubDir,
738  *this,
739  IOobject::NO_READ,
740  IOobject::AUTO_WRITE
741  ),
742  0
743  ),
744  clearedPrimitives_(false),
745  boundary_
746  (
747  IOobject
748  (
749  "boundary",
750  instance(),
751  meshSubDir,
752  *this,
753  IOobject::NO_READ,
754  IOobject::AUTO_WRITE
755  ),
756  *this,
757  boundaryFaces.size() + 1 // Add room for a default patch
758  ),
759  bounds_(points_, syncPar),
760  comm_(UPstream::worldComm),
761  geometricD_(Vector<label>::zero),
762  solutionD_(Vector<label>::zero),
763  tetBasePtIsPtr_(NULL),
764  cellTreePtr_(NULL),
765  pointZones_
766  (
767  IOobject
768  (
769  "pointZones",
770  instance(),
771  meshSubDir,
772  *this,
773  IOobject::NO_READ,
774  IOobject::NO_WRITE
775  ),
776  *this,
777  0
778  ),
779  faceZones_
780  (
781  IOobject
782  (
783  "faceZones",
784  instance(),
785  meshSubDir,
786  *this,
787  IOobject::NO_READ,
788  IOobject::NO_WRITE
789  ),
790  *this,
791  0
792  ),
793  cellZones_
794  (
795  IOobject
796  (
797  "cellZones",
798  instance(),
799  meshSubDir,
800  *this,
801  IOobject::NO_READ,
802  IOobject::NO_WRITE
803  ),
804  *this,
805  0
806  ),
807  globalMeshDataPtr_(NULL),
808  moving_(false),
809  topoChanging_(false),
810  curMotionTimeIndex_(time().timeIndex()),
811  oldPointsPtr_(NULL)
812 {
813  if (debug)
814  {
815  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
816  }
817 
818  // Remove all of the old mesh files if they exist
819  removeFiles(instance());
820 
821  // Calculate faces and cells
822  labelList patchSizes;
823  labelList patchStarts;
824  label defaultPatchStart;
825  label nFaces;
826  cellList cells;
827  setTopology
828  (
829  cellsAsShapes,
830  boundaryFaces,
831  boundaryPatchNames,
832  patchSizes,
833  patchStarts,
834  defaultPatchStart,
835  nFaces,
836  cells
837  );
838 
839  // Warning: Patches can only be added once the face list is
840  // completed, as they hold a subList of the face list
841  forAll(boundaryDicts, patchI)
842  {
843  dictionary patchDict(boundaryDicts[patchI]);
844 
845  patchDict.set("nFaces", patchSizes[patchI]);
846  patchDict.set("startFace", patchStarts[patchI]);
847 
848  // Add the patch to the list
849  boundary_.set
850  (
851  patchI,
853  (
854  boundaryPatchNames[patchI],
855  patchDict,
856  patchI,
857  boundary_
858  )
859  );
860  }
861 
862  label nAllPatches = boundaryFaces.size();
863 
864  label nDefaultFaces = nFaces - defaultPatchStart;
865  if (syncPar)
866  {
867  reduce(nDefaultFaces, sumOp<label>());
868  }
869 
870  if (nDefaultFaces > 0)
871  {
873  << "Found " << nDefaultFaces
874  << " undefined faces in mesh; adding to default patch." << endl;
875 
876  // Check if there already exists a defaultFaces patch as last patch
877  // and reuse it.
878  label patchI = findIndex(boundaryPatchNames, defaultBoundaryPatchName);
879 
880  if (patchI != -1)
881  {
882  if (patchI != boundaryFaces.size()-1 || boundary_[patchI].size())
883  {
885  << "Default patch " << boundary_[patchI].name()
886  << " already has faces in it or is not"
887  << " last in list of patches." << exit(FatalError);
888  }
889 
891  << "Reusing existing patch " << patchI
892  << " for undefined faces." << endl;
893 
894  boundary_.set
895  (
896  patchI,
898  (
899  boundary_[patchI].type(),
900  boundary_[patchI].name(),
901  nFaces - defaultPatchStart,
902  defaultPatchStart,
903  patchI,
904  boundary_
905  )
906  );
907  }
908  else
909  {
910  boundary_.set
911  (
912  nAllPatches,
914  (
915  defaultBoundaryPatchType,
916  defaultBoundaryPatchName,
917  nFaces - defaultPatchStart,
918  defaultPatchStart,
919  boundary_.size() - 1,
920  boundary_
921  )
922  );
923 
924  nAllPatches++;
925  }
926  }
927 
928  // Reset the size of the boundary
929  boundary_.setSize(nAllPatches);
930 
931  // Set the primitive mesh
932  initMesh(cells);
933 
934  if (syncPar)
935  {
936  // Calculate topology for the patches (processor-processor comms etc.)
937  boundary_.updateMesh();
938 
939  // Calculate the geometry for the patches (transformation tensors etc.)
940  boundary_.calcGeometry();
941  }
942 
943  if (debug)
944  {
945  if (checkMesh())
946  {
947  Info << "Mesh OK" << endl;
948  }
949  }
950 }
951 
952 
953 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::polyMeshGenChecks::checkMesh
bool checkMesh(const polyMeshGen &mesh, const bool report)
Check mesh for correctness. Returns false for no error.
Definition: polyMeshGenChecks.C:104
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
globalMeshData.H
Foam::polyMesh::facePatchFaceCells
labelList facePatchFaceCells(const faceList &patchFaces, const labelListList &pointCells, const faceListList &cellsFaceShapes, const label patchID) const
Definition: polyMeshFromShapeMesh.C:74
Foam::polyMesh::cellShapePointCells
labelListList cellShapePointCells(const cellShapeList &) const
Definition: polyMeshFromShapeMesh.C:37
indexedOctree.H
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::polyMesh::setTopology
void setTopology(const cellShapeList &cellsAsShapes, const faceListList &boundaryFaces, const wordList &boundaryPatchNames, labelList &patchSizes, labelList &patchStarts, label &defaultPatchStart, label &nFaces, cellList &cells)
Set faces_, calculate cells and patchStarts.
Definition: polyMeshFromShapeMesh.C:132
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
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::Info
messageStream Info
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
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
treeDataCell.H
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:57
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
DynamicList.H
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::polyMesh::polyMesh
polyMesh(const polyMesh &)
Disallow construct as copy.
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79