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, 2020 OpenFOAM Foundation
9  Copyright (C) 2018-2021 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  << " vertices " << UIndirectList<point>(points(), curFace)
122  << " does not have neighbour cell"
123  << " face: " << patchFaces[fI]
124  << abort(FatalError);
125  }
126  }
127 
128  return FaceCells;
129 }
130 
131 
132 void Foam::polyMesh::setTopology
133 (
134  const cellShapeList& cellsAsShapes,
135  const faceListList& boundaryFaces,
136  const wordList& boundaryPatchNames,
137  labelList& patchSizes,
138  labelList& patchStarts,
139  label& defaultPatchStart,
140  label& nFaces,
141  cellList& cells
142 )
143 {
144  // Calculate the faces of all cells
145  // Initialise maximum possible number of mesh faces to 0
146  label maxFaces = 0;
147 
148  // Set up a list of face shapes for each cell
149  faceListList cellsFaceShapes(cellsAsShapes.size());
150  cells.setSize(cellsAsShapes.size());
151 
152  forAll(cellsFaceShapes, celli)
153  {
154  cellsFaceShapes[celli] = cellsAsShapes[celli].faces();
155 
156  cells[celli].setSize(cellsFaceShapes[celli].size());
157 
158  // Initialise cells to -1 to flag undefined faces
159  static_cast<labelList&>(cells[celli]) = -1;
160 
161  // Count maximum possible number of mesh faces
162  maxFaces += cellsFaceShapes[celli].size();
163  }
164 
165  // Set size of faces array to maximum possible number of mesh faces
166  faces_.setSize(maxFaces);
167 
168  // Initialise number of faces to 0
169  nFaces = 0;
170 
171  // Set reference to point-cell addressing
172  labelListList PointCells = cellShapePointCells(cellsAsShapes);
173 
174  bool found = false;
175 
176  forAll(cells, celli)
177  {
178  // Note:
179  // Insertion cannot be done in one go as the faces need to be
180  // added into the list in the increasing order of neighbour
181  // cells. Therefore, all neighbours will be detected first
182  // and then added in the correct order.
183 
184  const faceList& curFaces = cellsFaceShapes[celli];
185 
186  // Record the neighbour cell
187  labelList neiCells(curFaces.size(), -1);
188 
189  // Record the face of neighbour cell
190  labelList faceOfNeiCell(curFaces.size(), -1);
191 
192  label nNeighbours = 0;
193 
194  // For all faces ...
195  forAll(curFaces, facei)
196  {
197  // Skip faces that have already been matched
198  if (cells[celli][facei] >= 0) continue;
199 
200  found = false;
201 
202  const face& curFace = curFaces[facei];
203 
204  // Get the list of labels
205  const labelList& curPoints = curFace;
206 
207  // For all points
208  forAll(curPoints, pointi)
209  {
210  // Get the list of cells sharing this point
211  const labelList& curNeighbours =
212  PointCells[curPoints[pointi]];
213 
214  // For all neighbours
215  forAll(curNeighbours, neiI)
216  {
217  label curNei = curNeighbours[neiI];
218 
219  // Reject neighbours with the lower label
220  if (curNei > celli)
221  {
222  // Get the list of search faces
223  const faceList& searchFaces = cellsFaceShapes[curNei];
224 
225  forAll(searchFaces, neiFacei)
226  {
227  if (searchFaces[neiFacei] == curFace)
228  {
229  // Match!!
230  found = true;
231 
232  // Record the neighbour cell and face
233  neiCells[facei] = curNei;
234  faceOfNeiCell[facei] = neiFacei;
235  nNeighbours++;
236 
237  break;
238  }
239  }
240  if (found) break;
241  }
242  if (found) break;
243  }
244  if (found) break;
245  } // End of current points
246  } // End of current faces
247 
248  // Add the faces in the increasing order of neighbours
249  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
250  {
251  // Find the lowest neighbour which is still valid
252  label nextNei = -1;
253  label minNei = cells.size();
254 
255  forAll(neiCells, ncI)
256  {
257  if (neiCells[ncI] > -1 && neiCells[ncI] < minNei)
258  {
259  nextNei = ncI;
260  minNei = neiCells[ncI];
261  }
262  }
263 
264  if (nextNei > -1)
265  {
266  // Add the face to the list of faces
267  faces_[nFaces] = curFaces[nextNei];
268 
269  // Set cell-face and cell-neighbour-face to current face label
270  cells[celli][nextNei] = nFaces;
271  cells[neiCells[nextNei]][faceOfNeiCell[nextNei]] = nFaces;
272 
273  // Stop the neighbour from being used again
274  neiCells[nextNei] = -1;
275 
276  // Increment number of faces counter
277  nFaces++;
278  }
279  else
280  {
282  << "Error in internal face insertion"
283  << abort(FatalError);
284  }
285  }
286  }
287 
288  // Do boundary faces
289  const label nInternalFaces = nFaces;
290 
291  patchSizes.setSize(boundaryFaces.size(), -1);
292  patchStarts.setSize(boundaryFaces.size(), -1);
293 
294  forAll(boundaryFaces, patchi)
295  {
296  const faceList& patchFaces = boundaryFaces[patchi];
297 
298  labelList curPatchFaceCells =
299  facePatchFaceCells
300  (
301  patchFaces,
302  PointCells,
303  cellsFaceShapes,
304  patchi
305  );
306 
307  // Grab the start label
308  label curPatchStart = nFaces;
309 
310  // Suppress multiple warnings per patch
311  bool patchWarned = false;
312 
313  forAll(patchFaces, facei)
314  {
315  const face& curFace = patchFaces[facei];
316 
317  const label cellInside = curPatchFaceCells[facei];
318 
319  // Get faces of the cell inside
320  const faceList& facesOfCellInside = cellsFaceShapes[cellInside];
321 
322  bool found = false;
323 
324  forAll(facesOfCellInside, cellFacei)
325  {
326  if (face::sameVertices(facesOfCellInside[cellFacei], curFace))
327  {
328  found = true;
329 
330  const label meshFacei = cells[cellInside][cellFacei];
331 
332  if (meshFacei >= 0)
333  {
334  // Already have mesh face for this side of the
335  // cellshape. This can happen for duplicate faces.
336  // It might be
337  // an error or explicitly desired (e.g. duplicate
338  // baffles or acmi). We could have a special 7-faced
339  // hex shape instead so we can have additional patches
340  // but that would be unworkable.
341  // So now either
342  // - exit with error
343  // - or warn and append face to addressing
344  // Note that duplicate baffles
345  // - cannot be on an internal faces
346  // - cannot be on the same patch (for now?)
347 
348  if
349  (
350  meshFacei < nInternalFaces
351  || meshFacei >= curPatchStart
352  )
353  {
355  << "Trying to specify a boundary face "
356  << curFace
357  << " on the face on cell " << cellInside
358  << " which is either an internal face"
359  << " or already belongs to the same patch."
360  << " This is face " << facei << " of patch "
361  << patchi << " named "
362  << boundaryPatchNames[patchi] << "."
363  << exit(FatalError);
364  }
365 
366 
367  if (!patchWarned)
368  {
370  << "Trying to specify a boundary face "
371  << curFace
372  << " on the face on cell " << cellInside
373  << " which is either an internal face"
374  << " or already belongs to some other patch."
375  << " This is face " << facei << " of patch "
376  << patchi << " named "
377  << boundaryPatchNames[patchi] << "."
378  //<< abort(FatalError);
379  << endl;
380  patchWarned = true;
381  }
382 
383  faces_.setSize(faces_.size()+1);
384 
385  // Set the patch face to corresponding cell-face
386  faces_[nFaces] = facesOfCellInside[cellFacei];
387 
388  cells[cellInside].append(nFaces);
389  }
390  else
391  {
392  // Set the patch face to corresponding cell-face
393  faces_[nFaces] = facesOfCellInside[cellFacei];
394 
395  cells[cellInside][cellFacei] = nFaces;
396  }
397 
398  break;
399  }
400  }
401 
402  if (!found)
403  {
405  << "face " << facei << " of patch " << patchi
406  << " does not seem to belong to cell " << cellInside
407  << " which, according to the addressing, "
408  << "should be next to it."
409  << abort(FatalError);
410  }
411 
412  // Increment the counter of faces
413  nFaces++;
414  }
415 
416  patchSizes[patchi] = nFaces - curPatchStart;
417  patchStarts[patchi] = curPatchStart;
418  }
419 
420  // Grab "non-existing" faces and put them into a default patch
421 
422  defaultPatchStart = nFaces;
423 
424  forAll(cells, celli)
425  {
426  labelList& curCellFaces = cells[celli];
427 
428  forAll(curCellFaces, facei)
429  {
430  if (curCellFaces[facei] == -1) // "non-existent" face
431  {
432  curCellFaces[facei] = nFaces;
433  faces_[nFaces] = cellsFaceShapes[celli][facei];
434 
435  nFaces++;
436  }
437  }
438  }
439 
440  // Reset the size of the face list
441  faces_.setSize(nFaces);
442 }
443 
444 
445 Foam::polyMesh::polyMesh
446 (
447  const IOobject& io,
448  pointField&& points,
449  const cellShapeList& cellsAsShapes,
450  const faceListList& boundaryFaces,
451  const wordList& boundaryPatchNames,
452  const wordList& boundaryPatchTypes,
453  const word& defaultBoundaryPatchName,
454  const word& defaultBoundaryPatchType,
455  const wordList& boundaryPatchPhysicalTypes,
456  const bool syncPar
457 )
458 :
459  objectRegistry(io),
460  primitiveMesh(),
461  points_
462  (
463  IOobject
464  (
465  "points",
466  instance(),
467  meshSubDir,
468  *this,
471  ),
472  std::move(points)
473  ),
474  faces_
475  (
476  IOobject
477  (
478  "faces",
479  instance(),
480  meshSubDir,
481  *this,
484  ),
485  0
486  ),
487  owner_
488  (
489  IOobject
490  (
491  "owner",
492  instance(),
493  meshSubDir,
494  *this,
497  ),
498  0
499  ),
500  neighbour_
501  (
502  IOobject
503  (
504  "neighbour",
505  instance(),
506  meshSubDir,
507  *this,
510  ),
511  0
512  ),
513  clearedPrimitives_(false),
514  boundary_
515  (
516  IOobject
517  (
518  "boundary",
519  instance(),
520  meshSubDir,
521  *this,
524  ),
525  *this,
526  boundaryFaces.size() + 1 // Add room for a default patch
527  ),
528  bounds_(points_, syncPar),
529  comm_(UPstream::worldComm),
530  geometricD_(Zero),
531  solutionD_(Zero),
532  tetBasePtIsPtr_(nullptr),
533  cellTreePtr_(nullptr),
534  pointZones_
535  (
536  IOobject
537  (
538  "pointZones",
539  instance(),
540  meshSubDir,
541  *this,
544  ),
545  *this,
546  0
547  ),
548  faceZones_
549  (
550  IOobject
551  (
552  "faceZones",
553  instance(),
554  meshSubDir,
555  *this,
558  ),
559  *this,
560  0
561  ),
562  cellZones_
563  (
564  IOobject
565  (
566  "cellZones",
567  instance(),
568  meshSubDir,
569  *this,
572  ),
573  *this,
574  0
575  ),
576  globalMeshDataPtr_(nullptr),
577  moving_(false),
578  topoChanging_(false),
579  storeOldCellCentres_(false),
580  curMotionTimeIndex_(time().timeIndex()),
581  oldPointsPtr_(nullptr),
582  oldCellCentresPtr_(nullptr)
583 {
584  DebugInfo
585  << "Constructing polyMesh from cell and boundary shapes." << endl;
586 
587  // Calculate faces and cells
588  labelList patchSizes;
589  labelList patchStarts;
590  label defaultPatchStart;
591  label nFaces;
592  cellList cells;
593  setTopology
594  (
595  cellsAsShapes,
596  boundaryFaces,
597  boundaryPatchNames,
598  patchSizes,
599  patchStarts,
600  defaultPatchStart,
601  nFaces,
602  cells
603  );
604 
605  // Warning: Patches can only be added once the face list is
606  // completed, as they hold a subList of the face list
607  forAll(boundaryFaces, patchi)
608  {
609  // Add the patch to the list
610  boundary_.set
611  (
612  patchi,
614  (
615  boundaryPatchTypes[patchi],
616  boundaryPatchNames[patchi],
617  patchSizes[patchi],
618  patchStarts[patchi],
619  patchi,
620  boundary_
621  )
622  );
623 
624  if
625  (
626  boundaryPatchPhysicalTypes.size()
627  && boundaryPatchPhysicalTypes[patchi].size()
628  )
629  {
630  boundary_[patchi].physicalType() =
631  boundaryPatchPhysicalTypes[patchi];
632  }
633  }
634 
635  label nAllPatches = boundaryFaces.size();
636 
637 
638  label nDefaultFaces = nFaces - defaultPatchStart;
639  if (syncPar)
640  {
641  reduce(nDefaultFaces, sumOp<label>());
642  }
643 
644  if (nDefaultFaces > 0)
645  {
647  << "Found " << nDefaultFaces
648  << " undefined faces in mesh; adding to default patch "
649  << defaultBoundaryPatchName << endl;
650 
651  // Check if there already exists a defaultFaces patch as last patch
652  // and reuse it.
653  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
654 
655  if (patchi != -1)
656  {
657  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
658  {
660  << "Default patch " << boundary_[patchi].name()
661  << " already has faces in it or is not"
662  << " last in list of patches." << exit(FatalError);
663  }
664 
666  << "Reusing existing patch " << patchi
667  << " for undefined faces." << endl;
668 
669  boundary_.set
670  (
671  patchi,
673  (
674  boundaryPatchTypes[patchi],
675  boundaryPatchNames[patchi],
676  nFaces - defaultPatchStart,
677  defaultPatchStart,
678  patchi,
679  boundary_
680  )
681  );
682  }
683  else
684  {
685  boundary_.set
686  (
687  nAllPatches,
689  (
690  defaultBoundaryPatchType,
691  defaultBoundaryPatchName,
692  nFaces - defaultPatchStart,
693  defaultPatchStart,
694  boundary_.size() - 1,
695  boundary_
696  )
697  );
698 
699  nAllPatches++;
700  }
701  }
702 
703  // Reset the size of the boundary
704  boundary_.setSize(nAllPatches);
705 
706  // Set the primitive mesh
707  initMesh(cells);
708 
709  if (syncPar)
710  {
711  // Calculate topology for the patches (processor-processor comms etc.)
712  boundary_.updateMesh();
713 
714  // Calculate the geometry for the patches (transformation tensors etc.)
715  boundary_.calcGeometry();
716  }
717 
718  if (debug)
719  {
720  if (checkMesh())
721  {
722  Info<< "Mesh OK" << endl;
723  }
724  }
725 }
726 
727 
728 Foam::polyMesh::polyMesh
729 (
730  const IOobject& io,
731  pointField&& points,
732  const cellShapeList& cellsAsShapes,
733  const faceListList& boundaryFaces,
734  const wordList& boundaryPatchNames,
736  const word& defaultBoundaryPatchName,
737  const word& defaultBoundaryPatchType,
738  const bool syncPar
739 )
740 :
741  objectRegistry(io),
742  primitiveMesh(),
743  points_
744  (
745  IOobject
746  (
747  "points",
748  instance(),
749  meshSubDir,
750  *this,
753  ),
754  std::move(points)
755  ),
756  faces_
757  (
758  IOobject
759  (
760  "faces",
761  instance(),
762  meshSubDir,
763  *this,
766  ),
767  0
768  ),
769  owner_
770  (
771  IOobject
772  (
773  "owner",
774  instance(),
775  meshSubDir,
776  *this,
779  ),
780  0
781  ),
782  neighbour_
783  (
784  IOobject
785  (
786  "neighbour",
787  instance(),
788  meshSubDir,
789  *this,
792  ),
793  0
794  ),
795  clearedPrimitives_(false),
796  boundary_
797  (
798  IOobject
799  (
800  "boundary",
801  instance(),
802  meshSubDir,
803  *this,
806  ),
807  *this,
808  boundaryFaces.size() + 1 // Add room for a default patch
809  ),
810  bounds_(points_, syncPar),
811  comm_(UPstream::worldComm),
812  geometricD_(Zero),
813  solutionD_(Zero),
814  tetBasePtIsPtr_(nullptr),
815  cellTreePtr_(nullptr),
816  pointZones_
817  (
818  IOobject
819  (
820  "pointZones",
821  instance(),
822  meshSubDir,
823  *this,
826  ),
827  *this,
828  0
829  ),
830  faceZones_
831  (
832  IOobject
833  (
834  "faceZones",
835  instance(),
836  meshSubDir,
837  *this,
840  ),
841  *this,
842  0
843  ),
844  cellZones_
845  (
846  IOobject
847  (
848  "cellZones",
849  instance(),
850  meshSubDir,
851  *this,
854  ),
855  *this,
856  0
857  ),
858  globalMeshDataPtr_(nullptr),
859  moving_(false),
860  topoChanging_(false),
861  storeOldCellCentres_(false),
862  curMotionTimeIndex_(time().timeIndex()),
863  oldPointsPtr_(nullptr),
864  oldCellCentresPtr_(nullptr)
865 {
866  DebugInfo
867  << "Constructing polyMesh from cell and boundary shapes." << endl;
868 
869  // Calculate faces and cells
870  labelList patchSizes;
871  labelList patchStarts;
872  label defaultPatchStart;
873  label nFaces;
874  cellList cells;
875  setTopology
876  (
877  cellsAsShapes,
878  boundaryFaces,
879  boundaryPatchNames,
880  patchSizes,
881  patchStarts,
882  defaultPatchStart,
883  nFaces,
884  cells
885  );
886 
887  // Warning: Patches can only be added once the face list is
888  // completed, as they hold a subList of the face list
889  forAll(boundaryDicts, patchi)
890  {
891  dictionary patchDict(boundaryDicts[patchi]);
892 
893  patchDict.set("nFaces", patchSizes[patchi]);
894  patchDict.set("startFace", patchStarts[patchi]);
895 
896  // Add the patch to the list
897  boundary_.set
898  (
899  patchi,
901  (
902  boundaryPatchNames[patchi],
903  patchDict,
904  patchi,
905  boundary_
906  )
907  );
908  }
909 
910  label nAllPatches = boundaryFaces.size();
911 
912  label nDefaultFaces = nFaces - defaultPatchStart;
913  if (syncPar)
914  {
915  reduce(nDefaultFaces, sumOp<label>());
916  }
917 
918  if (nDefaultFaces > 0)
919  {
921  << "Found " << nDefaultFaces
922  << " undefined faces in mesh; adding to default patch "
923  << defaultBoundaryPatchName << endl;
924 
925  // Check if there already exists a defaultFaces patch as last patch
926  // and reuse it.
927  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
928 
929  if (patchi != -1)
930  {
931  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
932  {
934  << "Default patch " << boundary_[patchi].name()
935  << " already has faces in it or is not"
936  << " last in list of patches." << exit(FatalError);
937  }
938 
940  << "Reusing existing patch " << patchi
941  << " for undefined faces." << endl;
942 
943  boundary_.set
944  (
945  patchi,
947  (
948  boundary_[patchi].type(),
949  boundary_[patchi].name(),
950  nFaces - defaultPatchStart,
951  defaultPatchStart,
952  patchi,
953  boundary_
954  )
955  );
956  }
957  else
958  {
959  boundary_.set
960  (
961  nAllPatches,
963  (
964  defaultBoundaryPatchType,
965  defaultBoundaryPatchName,
966  nFaces - defaultPatchStart,
967  defaultPatchStart,
968  boundary_.size() - 1,
969  boundary_
970  )
971  );
972 
973  nAllPatches++;
974  }
975  }
976 
977  // Reset the size of the boundary
978  boundary_.setSize(nAllPatches);
979 
980  // Set the primitive mesh
981  initMesh(cells);
982 
983  if (syncPar)
984  {
985  // Calculate topology for the patches (processor-processor comms etc.)
986  boundary_.updateMesh();
987 
988  // Calculate the geometry for the patches (transformation tensors etc.)
989  boundary_.calcGeometry();
990  }
991 
992  if (debug)
993  {
994  if (checkMesh())
995  {
996  Info<< "Mesh OK" << endl;
997  }
998  }
999 }
1000 
1001 
1002 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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:1069
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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:382
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
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:62
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
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
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:59
boundaryDicts
Info<< "Creating cells"<< endl;cellShapes=b.shapes();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:56
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:123
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 pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:35
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
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:453
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:382
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:293
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
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.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
DynamicList.H
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::IOobject::NO_READ
Definition: IOobject.H:188
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78