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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "polyMesh.H"
30 #include "Time.H"
31 #include "primitiveMesh.H"
32 #include "DynamicList.H"
33 #include "indexedOctree.H"
34 #include "treeDataCell.H"
35 #include "globalMeshData.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 Foam::labelListList Foam::polyMesh::cellShapePointCells
40 (
41  const cellShapeList& c
42 ) const
43 {
44  List<DynamicList<label>> pc(points().size());
45 
46  // For each cell
47  forAll(c, i)
48  {
49  // For each vertex
50  const labelList& labels = c[i];
51 
52  forAll(labels, j)
53  {
54  // Set working point label
55  label curPoint = labels[j];
56  DynamicList<label>& curPointCells = pc[curPoint];
57 
58  // Enter the cell label in the point's cell list
59  curPointCells.append(i);
60  }
61  }
62 
63  labelListList pointCellAddr(pc.size());
64 
65  forAll(pc, pointi)
66  {
67  pointCellAddr[pointi].transfer(pc[pointi]);
68  }
69 
70  return pointCellAddr;
71 }
72 
73 
74 Foam::labelList Foam::polyMesh::facePatchFaceCells
75 (
76  const faceList& patchFaces,
77  const labelListList& pointCells,
78  const faceListList& cellsFaceShapes,
79  const label patchID
80 ) const
81 {
82  bool found;
83 
84  labelList FaceCells(patchFaces.size());
85 
86  forAll(patchFaces, fI)
87  {
88  found = false;
89 
90  const face& curFace = patchFaces[fI];
91  const labelList& facePoints = patchFaces[fI];
92 
93  forAll(facePoints, pointi)
94  {
95  const labelList& facePointCells = pointCells[facePoints[pointi]];
96 
97  forAll(facePointCells, celli)
98  {
99  faceList cellFaces = cellsFaceShapes[facePointCells[celli]];
100 
101  forAll(cellFaces, cellFace)
102  {
103  if (face::sameVertices(cellFaces[cellFace], curFace))
104  {
105  // Found the cell corresponding to this face
106  FaceCells[fI] = facePointCells[celli];
107 
108  found = true;
109  }
110  if (found) break;
111  }
112  if (found) break;
113  }
114  if (found) break;
115  }
116 
117  if (!found)
118  {
120  << "face " << fI << " in patch " << patchID
121  << " does not have neighbour cell"
122  << " face: " << patchFaces[fI]
123  << abort(FatalError);
124  }
125  }
126 
127  return FaceCells;
128 }
129 
130 
131 void Foam::polyMesh::setTopology
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 number 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 number 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  // Get 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  const label nInternalFaces = nFaces;
289 
290  patchSizes.setSize(boundaryFaces.size(), -1);
291  patchStarts.setSize(boundaryFaces.size(), -1);
292 
293  forAll(boundaryFaces, patchi)
294  {
295  const faceList& patchFaces = boundaryFaces[patchi];
296 
297  labelList curPatchFaceCells =
298  facePatchFaceCells
299  (
300  patchFaces,
301  PointCells,
302  cellsFaceShapes,
303  patchi
304  );
305 
306  // Grab the start label
307  label curPatchStart = nFaces;
308 
309  // Suppress multiple warnings per patch
310  bool patchWarned = false;
311 
312  forAll(patchFaces, facei)
313  {
314  const face& curFace = patchFaces[facei];
315 
316  const label cellInside = curPatchFaceCells[facei];
317 
318  // Get faces of the cell inside
319  const faceList& facesOfCellInside = cellsFaceShapes[cellInside];
320 
321  bool found = false;
322 
323  forAll(facesOfCellInside, cellFacei)
324  {
325  if (face::sameVertices(facesOfCellInside[cellFacei], curFace))
326  {
327  found = true;
328 
329  const label meshFacei = cells[cellInside][cellFacei];
330 
331  if (meshFacei >= 0)
332  {
333  // Already have mesh face for this side of the
334  // cellshape. This can happen for duplicate faces.
335  // It might be
336  // an error or explicitly desired (e.g. duplicate
337  // baffles or acmi). We could have a special 7-faced
338  // hex shape instead so we can have additional patches
339  // but that would be unworkable.
340  // So now either
341  // - exit with error
342  // - or warn and append face to addressing
343  // Note that duplicate baffles
344  // - cannot be on an internal faces
345  // - cannot be on the same patch (for now?)
346 
347  if
348  (
349  meshFacei < nInternalFaces
350  || meshFacei >= curPatchStart
351  )
352  {
354  << "Trying to specify a boundary face "
355  << curFace
356  << " on the face on cell " << cellInside
357  << " which is either an internal face"
358  << " or already belongs to the same patch."
359  << " This is face " << facei << " of patch "
360  << patchi << " named "
361  << boundaryPatchNames[patchi] << "."
362  << exit(FatalError);
363  }
364 
365 
366  if (!patchWarned)
367  {
369  << "Trying to specify a boundary face "
370  << curFace
371  << " on the face on cell " << cellInside
372  << " which is either an internal face"
373  << " or already belongs to some other patch."
374  << " This is face " << facei << " of patch "
375  << patchi << " named "
376  << boundaryPatchNames[patchi] << "."
377  //<< abort(FatalError);
378  << endl;
379  patchWarned = true;
380  }
381 
382  faces_.setSize(faces_.size()+1);
383 
384  // Set the patch face to corresponding cell-face
385  faces_[nFaces] = facesOfCellInside[cellFacei];
386 
387  cells[cellInside].append(nFaces);
388  }
389  else
390  {
391  // Set the patch face to corresponding cell-face
392  faces_[nFaces] = facesOfCellInside[cellFacei];
393 
394  cells[cellInside][cellFacei] = nFaces;
395  }
396 
397  break;
398  }
399  }
400 
401  if (!found)
402  {
404  << "face " << facei << " of patch " << patchi
405  << " does not seem to belong to cell " << cellInside
406  << " which, according to the addressing, "
407  << "should be next to it."
408  << abort(FatalError);
409  }
410 
411  // Increment the counter of faces
412  nFaces++;
413  }
414 
415  patchSizes[patchi] = nFaces - curPatchStart;
416  patchStarts[patchi] = curPatchStart;
417  }
418 
419  // Grab "non-existing" faces and put them into a default patch
420 
421  defaultPatchStart = nFaces;
422 
423  forAll(cells, celli)
424  {
425  labelList& curCellFaces = cells[celli];
426 
427  forAll(curCellFaces, facei)
428  {
429  if (curCellFaces[facei] == -1) // "non-existent" face
430  {
431  curCellFaces[facei] = nFaces;
432  faces_[nFaces] = cellsFaceShapes[celli][facei];
433 
434  nFaces++;
435  }
436  }
437  }
438 
439  // Reset the size of the face list
440  faces_.setSize(nFaces);
441 }
442 
443 
444 Foam::polyMesh::polyMesh
445 (
446  const IOobject& io,
447  pointField&& points,
448  const cellShapeList& cellsAsShapes,
449  const faceListList& boundaryFaces,
450  const wordList& boundaryPatchNames,
451  const wordList& boundaryPatchTypes,
452  const word& defaultBoundaryPatchName,
453  const word& defaultBoundaryPatchType,
454  const wordList& boundaryPatchPhysicalTypes,
455  const bool syncPar
456 )
457 :
458  objectRegistry(io),
459  primitiveMesh(),
460  points_
461  (
462  IOobject
463  (
464  "points",
465  instance(),
466  meshSubDir,
467  *this,
470  ),
471  std::move(points)
472  ),
473  faces_
474  (
475  IOobject
476  (
477  "faces",
478  instance(),
479  meshSubDir,
480  *this,
483  ),
484  0
485  ),
486  owner_
487  (
488  IOobject
489  (
490  "owner",
491  instance(),
492  meshSubDir,
493  *this,
496  ),
497  0
498  ),
499  neighbour_
500  (
501  IOobject
502  (
503  "neighbour",
504  instance(),
505  meshSubDir,
506  *this,
509  ),
510  0
511  ),
512  clearedPrimitives_(false),
513  boundary_
514  (
515  IOobject
516  (
517  "boundary",
518  instance(),
519  meshSubDir,
520  *this,
523  ),
524  *this,
525  boundaryFaces.size() + 1 // Add room for a default patch
526  ),
527  bounds_(points_, syncPar),
528  comm_(UPstream::worldComm),
529  geometricD_(Zero),
530  solutionD_(Zero),
531  tetBasePtIsPtr_(nullptr),
532  cellTreePtr_(nullptr),
533  pointZones_
534  (
535  IOobject
536  (
537  "pointZones",
538  instance(),
539  meshSubDir,
540  *this,
543  ),
544  *this,
545  0
546  ),
547  faceZones_
548  (
549  IOobject
550  (
551  "faceZones",
552  instance(),
553  meshSubDir,
554  *this,
557  ),
558  *this,
559  0
560  ),
561  cellZones_
562  (
563  IOobject
564  (
565  "cellZones",
566  instance(),
567  meshSubDir,
568  *this,
571  ),
572  *this,
573  0
574  ),
575  globalMeshDataPtr_(nullptr),
576  moving_(false),
577  topoChanging_(false),
578  curMotionTimeIndex_(time().timeIndex()),
579  oldPointsPtr_(nullptr)
580 {
581  DebugInfo
582  << "Constructing polyMesh from cell and boundary shapes." << endl;
583 
584  // Calculate faces and cells
585  labelList patchSizes;
586  labelList patchStarts;
587  label defaultPatchStart;
588  label nFaces;
589  cellList cells;
590  setTopology
591  (
592  cellsAsShapes,
593  boundaryFaces,
594  boundaryPatchNames,
595  patchSizes,
596  patchStarts,
597  defaultPatchStart,
598  nFaces,
599  cells
600  );
601 
602  // Warning: Patches can only be added once the face list is
603  // completed, as they hold a subList of the face list
604  forAll(boundaryFaces, patchi)
605  {
606  // Add the patch to the list
607  boundary_.set
608  (
609  patchi,
611  (
612  boundaryPatchTypes[patchi],
613  boundaryPatchNames[patchi],
614  patchSizes[patchi],
615  patchStarts[patchi],
616  patchi,
617  boundary_
618  )
619  );
620 
621  if
622  (
623  boundaryPatchPhysicalTypes.size()
624  && boundaryPatchPhysicalTypes[patchi].size()
625  )
626  {
627  boundary_[patchi].physicalType() =
628  boundaryPatchPhysicalTypes[patchi];
629  }
630  }
631 
632  label nAllPatches = boundaryFaces.size();
633 
634 
635  label nDefaultFaces = nFaces - defaultPatchStart;
636  if (syncPar)
637  {
638  reduce(nDefaultFaces, sumOp<label>());
639  }
640 
641  if (nDefaultFaces > 0)
642  {
644  << "Found " << nDefaultFaces
645  << " undefined faces in mesh; adding to default patch "
646  << defaultBoundaryPatchName << endl;
647 
648  // Check if there already exists a defaultFaces patch as last patch
649  // and reuse it.
650  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
651 
652  if (patchi != -1)
653  {
654  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
655  {
657  << "Default patch " << boundary_[patchi].name()
658  << " already has faces in it or is not"
659  << " last in list of patches." << exit(FatalError);
660  }
661 
663  << "Reusing existing patch " << patchi
664  << " for undefined faces." << endl;
665 
666  boundary_.set
667  (
668  patchi,
670  (
671  boundaryPatchTypes[patchi],
672  boundaryPatchNames[patchi],
673  nFaces - defaultPatchStart,
674  defaultPatchStart,
675  patchi,
676  boundary_
677  )
678  );
679  }
680  else
681  {
682  boundary_.set
683  (
684  nAllPatches,
686  (
687  defaultBoundaryPatchType,
688  defaultBoundaryPatchName,
689  nFaces - defaultPatchStart,
690  defaultPatchStart,
691  boundary_.size() - 1,
692  boundary_
693  )
694  );
695 
696  nAllPatches++;
697  }
698  }
699 
700  // Reset the size of the boundary
701  boundary_.setSize(nAllPatches);
702 
703  // Set the primitive mesh
704  initMesh(cells);
705 
706  if (syncPar)
707  {
708  // Calculate topology for the patches (processor-processor comms etc.)
709  boundary_.updateMesh();
710 
711  // Calculate the geometry for the patches (transformation tensors etc.)
712  boundary_.calcGeometry();
713  }
714 
715  if (debug)
716  {
717  if (checkMesh())
718  {
719  Info<< "Mesh OK" << endl;
720  }
721  }
722 }
723 
724 
725 Foam::polyMesh::polyMesh
726 (
727  const IOobject& io,
728  pointField&& points,
729  const cellShapeList& cellsAsShapes,
730  const faceListList& boundaryFaces,
731  const wordList& boundaryPatchNames,
733  const word& defaultBoundaryPatchName,
734  const word& defaultBoundaryPatchType,
735  const bool syncPar
736 )
737 :
738  objectRegistry(io),
739  primitiveMesh(),
740  points_
741  (
742  IOobject
743  (
744  "points",
745  instance(),
746  meshSubDir,
747  *this,
750  ),
751  std::move(points)
752  ),
753  faces_
754  (
755  IOobject
756  (
757  "faces",
758  instance(),
759  meshSubDir,
760  *this,
763  ),
764  0
765  ),
766  owner_
767  (
768  IOobject
769  (
770  "owner",
771  instance(),
772  meshSubDir,
773  *this,
776  ),
777  0
778  ),
779  neighbour_
780  (
781  IOobject
782  (
783  "neighbour",
784  instance(),
785  meshSubDir,
786  *this,
789  ),
790  0
791  ),
792  clearedPrimitives_(false),
793  boundary_
794  (
795  IOobject
796  (
797  "boundary",
798  instance(),
799  meshSubDir,
800  *this,
803  ),
804  *this,
805  boundaryFaces.size() + 1 // Add room for a default patch
806  ),
807  bounds_(points_, syncPar),
808  comm_(UPstream::worldComm),
809  geometricD_(Zero),
810  solutionD_(Zero),
811  tetBasePtIsPtr_(nullptr),
812  cellTreePtr_(nullptr),
813  pointZones_
814  (
815  IOobject
816  (
817  "pointZones",
818  instance(),
819  meshSubDir,
820  *this,
823  ),
824  *this,
825  0
826  ),
827  faceZones_
828  (
829  IOobject
830  (
831  "faceZones",
832  instance(),
833  meshSubDir,
834  *this,
837  ),
838  *this,
839  0
840  ),
841  cellZones_
842  (
843  IOobject
844  (
845  "cellZones",
846  instance(),
847  meshSubDir,
848  *this,
851  ),
852  *this,
853  0
854  ),
855  globalMeshDataPtr_(nullptr),
856  moving_(false),
857  topoChanging_(false),
858  curMotionTimeIndex_(time().timeIndex()),
859  oldPointsPtr_(nullptr)
860 {
861  DebugInfo
862  << "Constructing polyMesh from cell and boundary shapes." << endl;
863 
864  // Calculate faces and cells
865  labelList patchSizes;
866  labelList patchStarts;
867  label defaultPatchStart;
868  label nFaces;
869  cellList cells;
870  setTopology
871  (
872  cellsAsShapes,
873  boundaryFaces,
874  boundaryPatchNames,
875  patchSizes,
876  patchStarts,
877  defaultPatchStart,
878  nFaces,
879  cells
880  );
881 
882  // Warning: Patches can only be added once the face list is
883  // completed, as they hold a subList of the face list
884  forAll(boundaryDicts, patchi)
885  {
886  dictionary patchDict(boundaryDicts[patchi]);
887 
888  patchDict.set("nFaces", patchSizes[patchi]);
889  patchDict.set("startFace", patchStarts[patchi]);
890 
891  // Add the patch to the list
892  boundary_.set
893  (
894  patchi,
896  (
897  boundaryPatchNames[patchi],
898  patchDict,
899  patchi,
900  boundary_
901  )
902  );
903  }
904 
905  label nAllPatches = boundaryFaces.size();
906 
907  label nDefaultFaces = nFaces - defaultPatchStart;
908  if (syncPar)
909  {
910  reduce(nDefaultFaces, sumOp<label>());
911  }
912 
913  if (nDefaultFaces > 0)
914  {
916  << "Found " << nDefaultFaces
917  << " undefined faces in mesh; adding to default patch "
918  << defaultBoundaryPatchName << endl;
919 
920  // Check if there already exists a defaultFaces patch as last patch
921  // and reuse it.
922  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
923 
924  if (patchi != -1)
925  {
926  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
927  {
929  << "Default patch " << boundary_[patchi].name()
930  << " already has faces in it or is not"
931  << " last in list of patches." << exit(FatalError);
932  }
933 
935  << "Reusing existing patch " << patchi
936  << " for undefined faces." << endl;
937 
938  boundary_.set
939  (
940  patchi,
942  (
943  boundary_[patchi].type(),
944  boundary_[patchi].name(),
945  nFaces - defaultPatchStart,
946  defaultPatchStart,
947  patchi,
948  boundary_
949  )
950  );
951  }
952  else
953  {
954  boundary_.set
955  (
956  nAllPatches,
958  (
959  defaultBoundaryPatchType,
960  defaultBoundaryPatchName,
961  nFaces - defaultPatchStart,
962  defaultPatchStart,
963  boundary_.size() - 1,
964  boundary_
965  )
966  );
967 
968  nAllPatches++;
969  }
970  }
971 
972  // Reset the size of the boundary
973  boundary_.setSize(nAllPatches);
974 
975  // Set the primitive mesh
976  initMesh(cells);
977 
978  if (syncPar)
979  {
980  // Calculate topology for the patches (processor-processor comms etc.)
981  boundary_.updateMesh();
982 
983  // Calculate the geometry for the patches (transformation tensors etc.)
984  boundary_.calcGeometry();
985  }
986 
987  if (debug)
988  {
989  if (checkMesh())
990  {
991  Info<< "Mesh OK" << endl;
992  }
993  }
994 }
995 
996 
997 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::HashTable< regIOobject * >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
globalMeshData.H
indexedOctree.H
Foam::cellShapeList
List< cellShape > cellShapeList
List of cellShapes and PtrList of List of cellShape.
Definition: cellShapeList.H:45
Foam::face::sameVertices
static bool sameVertices(const face &a, const face &b)
Return true if the faces have the same vertices.
Definition: face.C:402
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:847
polyMesh.H
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
boundaryDicts
Info<< "Creating boundary faces"<< endl;boundary.setSize(b.boundaryPatches().size());forAll(boundary, patchi) { faceList faces(b.boundaryPatches()[patchi].size());forAll(faces, facei) { faces[facei]=face(b.boundaryPatches()[patchi][facei]);} boundary[patchi].transfer(faces);} points.transfer(const_cast< pointField & >b.points()));}Info<< "Creating patch dictionaries"<< endl;wordList patchNames(boundary.size());forAll(patchNames, patchi){ patchNames[patchi]=polyPatch::defaultName(patchi);}PtrList< dictionary > boundaryDicts(boundary.size())
Definition: createBlockMesh.H:64
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::polyPatch::New
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:35
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
treeDataCell.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
found
bool found
Definition: TABSMDCalcMethod2.H:32
Time.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::faceListList
List< faceList > faceListList
A List of faceList.
Definition: faceListFwd.H:49
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:354
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::UPstream::worldComm
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:285
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:102
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
DynamicList.H
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::IOobject::NO_READ
Definition: IOobject.H:123
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78