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-2023 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& shapes
42 ) const
43 {
44  List<DynamicList<label>> pc(points().size());
45 
46  // For each cell
47  forAll(shapes, celli)
48  {
49  // For each cell vertex
50  for (const label pointi : shapes[celli])
51  {
52  // Enter the cell label in the point's cell list
53  pc[pointi].push_back(celli);
54  }
55  }
56 
57  labelListList pointCellAddr(pc.size());
58 
59  forAll(pc, pointi)
60  {
61  pointCellAddr[pointi].transfer(pc[pointi]);
62  }
63 
64  return pointCellAddr;
65 }
66 
67 
68 Foam::labelList Foam::polyMesh::facePatchFaceCells
69 (
70  const faceList& patchFaces,
71  const labelListList& pointCells,
72  const faceListList& cellsFaceShapes,
73  const label patchID
74 ) const
75 {
76  bool found;
77 
78  labelList FaceCells(patchFaces.size());
79 
80  forAll(patchFaces, fI)
81  {
82  found = false;
83 
84  const face& curFace = patchFaces[fI];
85  const labelList& facePoints = patchFaces[fI];
86 
87  forAll(facePoints, pointi)
88  {
89  const labelList& facePointCells = pointCells[facePoints[pointi]];
90 
91  forAll(facePointCells, celli)
92  {
93  faceList cellFaces = cellsFaceShapes[facePointCells[celli]];
94 
95  forAll(cellFaces, cellFace)
96  {
97  if (face::sameVertices(cellFaces[cellFace], curFace))
98  {
99  // Found the cell corresponding to this face
100  FaceCells[fI] = facePointCells[celli];
101 
102  found = true;
103  }
104  if (found) break;
105  }
106  if (found) break;
107  }
108  if (found) break;
109  }
110 
111  if (!found)
112  {
114  << "face " << fI << " in patch " << patchID
115  << " vertices " << UIndirectList<point>(points(), curFace)
116  << " does not have neighbour cell"
117  << " face: " << patchFaces[fI]
118  << abort(FatalError);
119  }
120  }
121 
122  return FaceCells;
123 }
124 
125 
126 void Foam::polyMesh::setTopology
127 (
128  const cellShapeList& cellsAsShapes,
129  const faceListList& boundaryFaces,
130  const wordList& boundaryPatchNames,
131  labelList& patchSizes,
132  labelList& patchStarts,
133  label& defaultPatchStart,
134  label& nFaces,
135  cellList& cells
136 )
137 {
138  // Calculate the faces of all cells
139  // Initialise maximum possible number of mesh faces to 0
140  label maxFaces = 0;
141 
142  // Set up a list of face shapes for each cell
143  faceListList cellsFaceShapes(cellsAsShapes.size());
144  cells.setSize(cellsAsShapes.size());
145 
146  forAll(cellsFaceShapes, celli)
147  {
148  cellsFaceShapes[celli] = cellsAsShapes[celli].faces();
149 
150  cells[celli].setSize(cellsFaceShapes[celli].size());
151 
152  // Initialise cells to -1 to flag undefined faces
153  static_cast<labelList&>(cells[celli]) = -1;
154 
155  // Count maximum possible number of mesh faces
156  maxFaces += cellsFaceShapes[celli].size();
157  }
158 
159  // Set size of faces array to maximum possible number of mesh faces
160  faces_.setSize(maxFaces);
161 
162  // Initialise number of faces to 0
163  nFaces = 0;
164 
165  // Set reference to point-cell addressing
166  labelListList PointCells = cellShapePointCells(cellsAsShapes);
167 
168  bool found = false;
169 
170  forAll(cells, celli)
171  {
172  // Note:
173  // Insertion cannot be done in one go as the faces need to be
174  // added into the list in the increasing order of neighbour
175  // cells. Therefore, all neighbours will be detected first
176  // and then added in the correct order.
177 
178  const faceList& curFaces = cellsFaceShapes[celli];
179 
180  // Record the neighbour cell
181  labelList neiCells(curFaces.size(), -1);
182 
183  // Record the face of neighbour cell
184  labelList faceOfNeiCell(curFaces.size(), -1);
185 
186  label nNeighbours = 0;
187 
188  // For all faces ...
189  forAll(curFaces, facei)
190  {
191  // Skip faces that have already been matched
192  if (cells[celli][facei] >= 0) continue;
193 
194  found = false;
195 
196  const face& curFace = curFaces[facei];
197 
198  // Get the list of labels
199  const labelList& curPoints = curFace;
200 
201  // For all points
202  forAll(curPoints, pointi)
203  {
204  // Get the list of cells sharing this point
205  const labelList& curNeighbours =
206  PointCells[curPoints[pointi]];
207 
208  // For all neighbours
209  forAll(curNeighbours, neiI)
210  {
211  label curNei = curNeighbours[neiI];
212 
213  // Reject neighbours with the lower label
214  if (curNei > celli)
215  {
216  // Get the list of search faces
217  const faceList& searchFaces = cellsFaceShapes[curNei];
218 
219  forAll(searchFaces, neiFacei)
220  {
221  if (searchFaces[neiFacei] == curFace)
222  {
223  // Match!!
224  found = true;
225 
226  // Record the neighbour cell and face
227  neiCells[facei] = curNei;
228  faceOfNeiCell[facei] = neiFacei;
229  nNeighbours++;
230 
231  break;
232  }
233  }
234  if (found) break;
235  }
236  if (found) break;
237  }
238  if (found) break;
239  } // End of current points
240  } // End of current faces
241 
242  // Add the faces in the increasing order of neighbours
243  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
244  {
245  // Find the lowest neighbour which is still valid
246  label nextNei = -1;
247  label minNei = cells.size();
248 
249  forAll(neiCells, ncI)
250  {
251  if (neiCells[ncI] > -1 && neiCells[ncI] < minNei)
252  {
253  nextNei = ncI;
254  minNei = neiCells[ncI];
255  }
256  }
257 
258  if (nextNei > -1)
259  {
260  // Add the face to the list of faces
261  faces_[nFaces] = curFaces[nextNei];
262 
263  // Set cell-face and cell-neighbour-face to current face label
264  cells[celli][nextNei] = nFaces;
265  cells[neiCells[nextNei]][faceOfNeiCell[nextNei]] = nFaces;
266 
267  // Stop the neighbour from being used again
268  neiCells[nextNei] = -1;
269 
270  // Increment number of faces counter
271  nFaces++;
272  }
273  else
274  {
276  << "Error in internal face insertion"
277  << abort(FatalError);
278  }
279  }
280  }
281 
282  // Do boundary faces
283  const label nInternalFaces = nFaces;
284 
285  patchSizes.setSize(boundaryFaces.size(), -1);
286  patchStarts.setSize(boundaryFaces.size(), -1);
287 
288  forAll(boundaryFaces, patchi)
289  {
290  const faceList& patchFaces = boundaryFaces[patchi];
291 
292  labelList curPatchFaceCells =
293  facePatchFaceCells
294  (
295  patchFaces,
296  PointCells,
297  cellsFaceShapes,
298  patchi
299  );
300 
301  // Grab the start label
302  label curPatchStart = nFaces;
303 
304  // Suppress multiple warnings per patch
305  bool patchWarned = false;
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  found = true;
323 
324  const label meshFacei = cells[cellInside][cellFacei];
325 
326  if (meshFacei >= 0)
327  {
328  // Already have mesh face for this side of the
329  // cellshape. This can happen for duplicate faces.
330  // It might be
331  // an error or explicitly desired (e.g. duplicate
332  // baffles or acmi). We could have a special 7-faced
333  // hex shape instead so we can have additional patches
334  // but that would be unworkable.
335  // So now either
336  // - exit with error
337  // - or warn and append face to addressing
338  // Note that duplicate baffles
339  // - cannot be on an internal faces
340  // - cannot be on the same patch (for now?)
341 
342  if
343  (
344  meshFacei < nInternalFaces
345  || meshFacei >= curPatchStart
346  )
347  {
349  << "Trying to specify a boundary face "
350  << curFace
351  << " on the face on cell " << cellInside
352  << " which is either an internal face"
353  << " or already belongs to the same patch."
354  << " This is face " << facei << " of patch "
355  << patchi << " named "
356  << boundaryPatchNames[patchi] << "."
357  << exit(FatalError);
358  }
359 
360 
361  if (!patchWarned)
362  {
364  << "Trying to specify a boundary face "
365  << curFace
366  << " on the face on cell " << cellInside
367  << " which is either an internal face"
368  << " or already belongs to some other patch."
369  << " This is face " << facei << " of patch "
370  << patchi << " named "
371  << boundaryPatchNames[patchi] << "."
372  //<< abort(FatalError);
373  << endl;
374  patchWarned = true;
375  }
376 
377  faces_.setSize(faces_.size()+1);
378 
379  // Set the patch face to corresponding cell-face
380  faces_[nFaces] = facesOfCellInside[cellFacei];
381 
382  cells[cellInside].append(nFaces);
383  }
384  else
385  {
386  // Set the patch face to corresponding cell-face
387  faces_[nFaces] = facesOfCellInside[cellFacei];
388 
389  cells[cellInside][cellFacei] = nFaces;
390  }
391 
392  break;
393  }
394  }
395 
396  if (!found)
397  {
399  << "face " << facei << " of patch " << patchi
400  << " does not seem to belong to cell " << cellInside
401  << " which, according to the addressing, "
402  << "should be next to it."
403  << abort(FatalError);
404  }
405 
406  // Increment the counter of faces
407  nFaces++;
408  }
409 
410  patchSizes[patchi] = nFaces - curPatchStart;
411  patchStarts[patchi] = curPatchStart;
412  }
413 
414  // Grab "non-existing" faces and put them into a default patch
415 
416  defaultPatchStart = nFaces;
417 
418  forAll(cells, celli)
419  {
420  labelList& curCellFaces = cells[celli];
421 
422  forAll(curCellFaces, facei)
423  {
424  if (curCellFaces[facei] == -1) // "non-existent" face
425  {
426  curCellFaces[facei] = nFaces;
427  faces_[nFaces] = cellsFaceShapes[celli][facei];
428 
429  nFaces++;
430  }
431  }
432  }
434  // Reset the size of the face list
435  faces_.setSize(nFaces);
436 }
437 
438 
439 Foam::polyMesh::polyMesh
440 (
441  const IOobject& io,
442  pointField&& points,
443  const cellShapeList& cellsAsShapes,
444  const faceListList& boundaryFaces,
445  const wordList& boundaryPatchNames,
446  const wordList& boundaryPatchTypes,
447  const word& defaultBoundaryPatchName,
448  const word& defaultBoundaryPatchType,
449  const wordList& boundaryPatchPhysicalTypes,
450  const bool syncPar
451 )
452 :
454  primitiveMesh(),
455  data_(static_cast<const objectRegistry&>(*this)),
456  points_
457  (
458  IOobject
459  (
460  "points",
461  instance(),
462  meshSubDir,
463  *this,
464  IOobject::NO_READ,
465  IOobject::AUTO_WRITE
466  ),
467  std::move(points)
468  ),
469  faces_
470  (
471  IOobject
472  (
473  "faces",
474  instance(),
475  meshSubDir,
476  *this,
477  IOobject::NO_READ,
478  IOobject::AUTO_WRITE
479  ),
480  Foam::zero{}
481  ),
482  owner_
483  (
484  IOobject
485  (
486  "owner",
487  instance(),
488  meshSubDir,
489  *this,
492  ),
493  Foam::zero{}
494  ),
495  neighbour_
496  (
497  IOobject
498  (
499  "neighbour",
500  instance(),
501  meshSubDir,
502  *this,
505  ),
506  Foam::zero{}
507  ),
508  clearedPrimitives_(false),
509  boundary_
510  (
511  IOobject
512  (
513  "boundary",
514  instance(),
515  meshSubDir,
516  *this,
519  ),
520  *this,
521  boundaryFaces.size() + 1 // Add room for a default patch
522  ),
523  bounds_(points_, syncPar),
524  comm_(UPstream::worldComm),
525  geometricD_(Zero),
526  solutionD_(Zero),
527  tetBasePtIsPtr_(nullptr),
528  cellTreePtr_(nullptr),
529  pointZones_
530  (
531  IOobject
532  (
533  "pointZones",
534  instance(),
535  meshSubDir,
536  *this,
539  ),
540  *this,
541  Foam::zero{}
542  ),
543  faceZones_
544  (
545  IOobject
546  (
547  "faceZones",
548  instance(),
549  meshSubDir,
550  *this,
553  ),
554  *this,
555  Foam::zero{}
556  ),
557  cellZones_
558  (
559  IOobject
560  (
561  "cellZones",
562  instance(),
563  meshSubDir,
564  *this,
567  ),
568  *this,
569  Foam::zero{}
570  ),
571  globalMeshDataPtr_(nullptr),
572  moving_(false),
573  topoChanging_(false),
574  storeOldCellCentres_(false),
575  curMotionTimeIndex_(time().timeIndex()),
576  oldPointsPtr_(nullptr),
577  oldCellCentresPtr_(nullptr)
578 {
579  DebugInfo
580  << "Constructing polyMesh from cell and boundary shapes." << endl;
581 
582  // Calculate faces and cells
583  labelList patchSizes;
584  labelList patchStarts;
585  label defaultPatchStart;
586  label nFaces;
587  cellList cells;
588  setTopology
589  (
590  cellsAsShapes,
591  boundaryFaces,
592  boundaryPatchNames,
593  patchSizes,
594  patchStarts,
595  defaultPatchStart,
596  nFaces,
597  cells
598  );
599 
600  // Warning: Patches can only be added once the face list is
601  // completed, as they hold a subList of the face list
602  forAll(boundaryFaces, patchi)
603  {
604  // Add the patch to the list
605  boundary_.set
606  (
607  patchi,
609  (
610  boundaryPatchTypes[patchi],
611  boundaryPatchNames[patchi],
612  patchSizes[patchi],
613  patchStarts[patchi],
614  patchi,
615  boundary_
616  )
617  );
618 
619  if
620  (
621  boundaryPatchPhysicalTypes.size()
622  && boundaryPatchPhysicalTypes[patchi].size()
623  )
624  {
625  boundary_[patchi].physicalType() =
626  boundaryPatchPhysicalTypes[patchi];
627  }
628  }
629 
630  label nAllPatches = boundaryFaces.size();
631 
632  label nDefaultFaces = nFaces - defaultPatchStart;
633  if (syncPar)
634  {
635  reduce(nDefaultFaces, sumOp<label>());
636  }
637 
638  if (nDefaultFaces > 0)
639  {
641  << "Found " << nDefaultFaces
642  << " undefined faces in mesh; adding to default patch "
643  << defaultBoundaryPatchName << endl;
644 
645  // Check if there already exists a defaultFaces patch as last patch
646  // and reuse it.
647  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
648 
649  if (patchi != -1)
650  {
651  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
652  {
654  << "Default patch " << boundary_[patchi].name()
655  << " already has faces in it or is not"
656  << " last in list of patches." << exit(FatalError);
657  }
658 
660  << "Reusing existing patch " << patchi
661  << " for undefined faces." << endl;
662 
663  boundary_.set
664  (
665  patchi,
667  (
668  boundaryPatchTypes[patchi],
669  boundaryPatchNames[patchi],
670  nFaces - defaultPatchStart,
671  defaultPatchStart,
672  patchi,
673  boundary_
674  )
675  );
676  }
677  else
678  {
679  boundary_.set
680  (
681  nAllPatches,
683  (
684  defaultBoundaryPatchType,
685  defaultBoundaryPatchName,
686  nFaces - defaultPatchStart,
687  defaultPatchStart,
688  boundary_.size() - 1,
689  boundary_
690  )
691  );
692 
693  nAllPatches++;
694  }
695  }
696 
697  // Reset the size of the boundary
698  boundary_.setSize(nAllPatches);
699 
700  // Set the primitive mesh
701  initMesh(cells);
702 
703  if (syncPar)
704  {
705  // Calculate topology for the patches (processor-processor comms etc.)
706  boundary_.updateMesh();
707 
708  // Calculate the geometry for the patches (transformation tensors etc.)
709  boundary_.calcGeometry();
710  }
711 
712  if (debug)
713  {
714  if (checkMesh())
715  {
716  Info<< "Mesh OK" << endl;
717  }
718  }
719 }
720 
721 
722 Foam::polyMesh::polyMesh
723 (
724  const IOobject& io,
725  pointField&& points,
726  const cellShapeList& cellsAsShapes,
727  const faceListList& boundaryFaces,
728  const wordList& boundaryPatchNames,
729  const PtrList<dictionary>& boundaryDicts,
730  const word& defaultBoundaryPatchName,
731  const word& defaultBoundaryPatchType,
732  const bool syncPar
733 )
734 :
735  objectRegistry(io),
736  primitiveMesh(),
737  data_(static_cast<const objectRegistry&>(*this)),
738  points_
739  (
740  IOobject
741  (
742  "points",
743  instance(),
744  meshSubDir,
745  *this,
746  IOobject::NO_READ,
747  IOobject::AUTO_WRITE
748  ),
749  std::move(points)
750  ),
751  faces_
752  (
753  IOobject
754  (
755  "faces",
756  instance(),
757  meshSubDir,
758  *this,
759  IOobject::NO_READ,
760  IOobject::AUTO_WRITE
761  ),
762  Foam::zero{}
763  ),
764  owner_
765  (
766  IOobject
767  (
768  "owner",
769  instance(),
770  meshSubDir,
771  *this,
774  ),
775  Foam::zero{}
776  ),
777  neighbour_
778  (
779  IOobject
780  (
781  "neighbour",
782  instance(),
783  meshSubDir,
784  *this,
787  ),
788  Foam::zero{}
789  ),
790  clearedPrimitives_(false),
791  boundary_
792  (
793  IOobject
794  (
795  "boundary",
796  instance(),
797  meshSubDir,
798  *this,
801  ),
802  *this,
803  boundaryFaces.size() + 1 // Add room for a default patch
804  ),
805  bounds_(points_, syncPar),
806  comm_(UPstream::worldComm),
807  geometricD_(Zero),
808  solutionD_(Zero),
809  tetBasePtIsPtr_(nullptr),
810  cellTreePtr_(nullptr),
811  pointZones_
812  (
813  IOobject
814  (
815  "pointZones",
816  instance(),
817  meshSubDir,
818  *this,
821  ),
822  *this,
823  Foam::zero{}
824  ),
825  faceZones_
826  (
827  IOobject
828  (
829  "faceZones",
830  instance(),
831  meshSubDir,
832  *this,
835  ),
836  *this,
837  Foam::zero{}
838  ),
839  cellZones_
840  (
841  IOobject
842  (
843  "cellZones",
844  instance(),
845  meshSubDir,
846  *this,
849  ),
850  *this,
851  Foam::zero{}
852  ),
853  globalMeshDataPtr_(nullptr),
854  moving_(false),
855  topoChanging_(false),
856  storeOldCellCentres_(false),
857  curMotionTimeIndex_(time().timeIndex()),
858  oldPointsPtr_(nullptr),
859  oldCellCentresPtr_(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 // ************************************************************************* //
List< faceList > faceListList
List of faceList.
Definition: faceListFwd.H:44
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
List< cellShape > cellShapeList
List of cellShape.
Definition: cellShapeList.H:32
Ignore writing from objectRegistry::writeObject()
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:409
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:316
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
const pointField & points
A class for handling words, derived from Foam::string.
Definition: word.H:63
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
int debug
Static debugging option.
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())
List< word > wordList
List of word.
Definition: fileName.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
Nothing to be read.
Automatically write from objectRegistry::writeObject()
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
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:28
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Registry of regIOobjects.
static bool sameVertices(const face &a, const face &b)
True if the faces have all the same vertices.
Definition: face.C:374
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127