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 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 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  // 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 
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 
392 Foam::polyMesh::polyMesh
393 (
394  const IOobject& io,
395  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,
418  ),
419  std::move(points)
420  ),
421  faces_
422  (
423  IOobject
424  (
425  "faces",
426  instance(),
427  meshSubDir,
428  *this,
431  ),
432  0
433  ),
434  owner_
435  (
436  IOobject
437  (
438  "owner",
439  instance(),
440  meshSubDir,
441  *this,
444  ),
445  0
446  ),
447  neighbour_
448  (
449  IOobject
450  (
451  "neighbour",
452  instance(),
453  meshSubDir,
454  *this,
457  ),
458  0
459  ),
460  clearedPrimitives_(false),
461  boundary_
462  (
463  IOobject
464  (
465  "boundary",
466  instance(),
467  meshSubDir,
468  *this,
471  ),
472  *this,
473  boundaryFaces.size() + 1 // Add room for a default patch
474  ),
475  bounds_(points_, syncPar),
476  comm_(UPstream::worldComm),
477  geometricD_(Zero),
478  solutionD_(Zero),
479  tetBasePtIsPtr_(nullptr),
480  cellTreePtr_(nullptr),
481  pointZones_
482  (
483  IOobject
484  (
485  "pointZones",
486  instance(),
487  meshSubDir,
488  *this,
491  ),
492  *this,
493  0
494  ),
495  faceZones_
496  (
497  IOobject
498  (
499  "faceZones",
500  instance(),
501  meshSubDir,
502  *this,
505  ),
506  *this,
507  0
508  ),
509  cellZones_
510  (
511  IOobject
512  (
513  "cellZones",
514  instance(),
515  meshSubDir,
516  *this,
519  ),
520  *this,
521  0
522  ),
523  globalMeshDataPtr_(nullptr),
524  moving_(false),
525  topoChanging_(false),
526  curMotionTimeIndex_(time().timeIndex()),
527  oldPointsPtr_(nullptr)
528 {
529  if (debug)
530  {
531  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
532  }
533 
534  // Calculate faces and cells
535  labelList patchSizes;
536  labelList patchStarts;
537  label defaultPatchStart;
538  label nFaces;
539  cellList cells;
540  setTopology
541  (
542  cellsAsShapes,
543  boundaryFaces,
544  boundaryPatchNames,
545  patchSizes,
546  patchStarts,
547  defaultPatchStart,
548  nFaces,
549  cells
550  );
551 
552  // Warning: Patches can only be added once the face list is
553  // completed, as they hold a subList of the face list
554  forAll(boundaryFaces, patchi)
555  {
556  // Add the patch to the list
557  boundary_.set
558  (
559  patchi,
561  (
562  boundaryPatchTypes[patchi],
563  boundaryPatchNames[patchi],
564  patchSizes[patchi],
565  patchStarts[patchi],
566  patchi,
567  boundary_
568  )
569  );
570 
571  if
572  (
573  boundaryPatchPhysicalTypes.size()
574  && boundaryPatchPhysicalTypes[patchi].size()
575  )
576  {
577  boundary_[patchi].physicalType() =
578  boundaryPatchPhysicalTypes[patchi];
579  }
580  }
581 
582  label nAllPatches = boundaryFaces.size();
583 
584 
585  label nDefaultFaces = nFaces - defaultPatchStart;
586  if (syncPar)
587  {
588  reduce(nDefaultFaces, sumOp<label>());
589  }
590 
591  if (nDefaultFaces > 0)
592  {
594  << "Found " << nDefaultFaces
595  << " undefined faces in mesh; adding to default patch "
596  << defaultBoundaryPatchName << endl;
597 
598  // Check if there already exists a defaultFaces patch as last patch
599  // and reuse it.
600  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
601 
602  if (patchi != -1)
603  {
604  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
605  {
607  << "Default patch " << boundary_[patchi].name()
608  << " already has faces in it or is not"
609  << " last in list of patches." << exit(FatalError);
610  }
611 
613  << "Reusing existing patch " << patchi
614  << " for undefined faces." << endl;
615 
616  boundary_.set
617  (
618  patchi,
620  (
621  boundaryPatchTypes[patchi],
622  boundaryPatchNames[patchi],
623  nFaces - defaultPatchStart,
624  defaultPatchStart,
625  patchi,
626  boundary_
627  )
628  );
629  }
630  else
631  {
632  boundary_.set
633  (
634  nAllPatches,
636  (
637  defaultBoundaryPatchType,
638  defaultBoundaryPatchName,
639  nFaces - defaultPatchStart,
640  defaultPatchStart,
641  boundary_.size() - 1,
642  boundary_
643  )
644  );
645 
646  nAllPatches++;
647  }
648  }
649 
650  // Reset the size of the boundary
651  boundary_.setSize(nAllPatches);
652 
653  // Set the primitive mesh
654  initMesh(cells);
655 
656  if (syncPar)
657  {
658  // Calculate topology for the patches (processor-processor comms etc.)
659  boundary_.updateMesh();
660 
661  // Calculate the geometry for the patches (transformation tensors etc.)
662  boundary_.calcGeometry();
663  }
664 
665  if (debug)
666  {
667  if (checkMesh())
668  {
669  Info<< "Mesh OK" << endl;
670  }
671  }
672 }
673 
674 
675 Foam::polyMesh::polyMesh
676 (
677  const IOobject& io,
678  pointField&& points,
679  const cellShapeList& cellsAsShapes,
680  const faceListList& boundaryFaces,
681  const wordList& boundaryPatchNames,
683  const word& defaultBoundaryPatchName,
684  const word& defaultBoundaryPatchType,
685  const bool syncPar
686 )
687 :
688  objectRegistry(io),
689  primitiveMesh(),
690  points_
691  (
692  IOobject
693  (
694  "points",
695  instance(),
696  meshSubDir,
697  *this,
700  ),
701  std::move(points)
702  ),
703  faces_
704  (
705  IOobject
706  (
707  "faces",
708  instance(),
709  meshSubDir,
710  *this,
713  ),
714  0
715  ),
716  owner_
717  (
718  IOobject
719  (
720  "owner",
721  instance(),
722  meshSubDir,
723  *this,
726  ),
727  0
728  ),
729  neighbour_
730  (
731  IOobject
732  (
733  "neighbour",
734  instance(),
735  meshSubDir,
736  *this,
739  ),
740  0
741  ),
742  clearedPrimitives_(false),
743  boundary_
744  (
745  IOobject
746  (
747  "boundary",
748  instance(),
749  meshSubDir,
750  *this,
753  ),
754  *this,
755  boundaryFaces.size() + 1 // Add room for a default patch
756  ),
757  bounds_(points_, syncPar),
758  comm_(UPstream::worldComm),
759  geometricD_(Zero),
760  solutionD_(Zero),
761  tetBasePtIsPtr_(nullptr),
762  cellTreePtr_(nullptr),
763  pointZones_
764  (
765  IOobject
766  (
767  "pointZones",
768  instance(),
769  meshSubDir,
770  *this,
773  ),
774  *this,
775  0
776  ),
777  faceZones_
778  (
779  IOobject
780  (
781  "faceZones",
782  instance(),
783  meshSubDir,
784  *this,
787  ),
788  *this,
789  0
790  ),
791  cellZones_
792  (
793  IOobject
794  (
795  "cellZones",
796  instance(),
797  meshSubDir,
798  *this,
801  ),
802  *this,
803  0
804  ),
805  globalMeshDataPtr_(nullptr),
806  moving_(false),
807  topoChanging_(false),
808  curMotionTimeIndex_(time().timeIndex()),
809  oldPointsPtr_(nullptr)
810 {
811  if (debug)
812  {
813  Info<<"Constructing polyMesh from cell and boundary shapes." << endl;
814  }
815 
816  // Calculate faces and cells
817  labelList patchSizes;
818  labelList patchStarts;
819  label defaultPatchStart;
820  label nFaces;
821  cellList cells;
822  setTopology
823  (
824  cellsAsShapes,
825  boundaryFaces,
826  boundaryPatchNames,
827  patchSizes,
828  patchStarts,
829  defaultPatchStart,
830  nFaces,
831  cells
832  );
833 
834  // Warning: Patches can only be added once the face list is
835  // completed, as they hold a subList of the face list
836  forAll(boundaryDicts, patchi)
837  {
838  dictionary patchDict(boundaryDicts[patchi]);
839 
840  patchDict.set("nFaces", patchSizes[patchi]);
841  patchDict.set("startFace", patchStarts[patchi]);
842 
843  // Add the patch to the list
844  boundary_.set
845  (
846  patchi,
848  (
849  boundaryPatchNames[patchi],
850  patchDict,
851  patchi,
852  boundary_
853  )
854  );
855  }
856 
857  label nAllPatches = boundaryFaces.size();
858 
859  label nDefaultFaces = nFaces - defaultPatchStart;
860  if (syncPar)
861  {
862  reduce(nDefaultFaces, sumOp<label>());
863  }
864 
865  if (nDefaultFaces > 0)
866  {
868  << "Found " << nDefaultFaces
869  << " undefined faces in mesh; adding to default patch "
870  << defaultBoundaryPatchName << endl;
871 
872  // Check if there already exists a defaultFaces patch as last patch
873  // and reuse it.
874  label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
875 
876  if (patchi != -1)
877  {
878  if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
879  {
881  << "Default patch " << boundary_[patchi].name()
882  << " already has faces in it or is not"
883  << " last in list of patches." << exit(FatalError);
884  }
885 
887  << "Reusing existing patch " << patchi
888  << " for undefined faces." << endl;
889 
890  boundary_.set
891  (
892  patchi,
894  (
895  boundary_[patchi].type(),
896  boundary_[patchi].name(),
897  nFaces - defaultPatchStart,
898  defaultPatchStart,
899  patchi,
900  boundary_
901  )
902  );
903  }
904  else
905  {
906  boundary_.set
907  (
908  nAllPatches,
910  (
911  defaultBoundaryPatchType,
912  defaultBoundaryPatchName,
913  nFaces - defaultPatchStart,
914  defaultPatchStart,
915  boundary_.size() - 1,
916  boundary_
917  )
918  );
919 
920  nAllPatches++;
921  }
922  }
923 
924  // Reset the size of the boundary
925  boundary_.setSize(nAllPatches);
926 
927  // Set the primitive mesh
928  initMesh(cells);
929 
930  if (syncPar)
931  {
932  // Calculate topology for the patches (processor-processor comms etc.)
933  boundary_.updateMesh();
934 
935  // Calculate the geometry for the patches (transformation tensors etc.)
936  boundary_.calcGeometry();
937  }
938 
939  if (debug)
940  {
941  if (checkMesh())
942  {
943  Info << "Mesh OK" << endl;
944  }
945  }
946 }
947 
948 
949 // ************************************************************************* //
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:74
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.
Definition: zero.H:128
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
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:842
polyMesh.H
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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
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:62
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:65
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
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:355
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]="patch"+Foam::name(patchi);}PtrList< dictionary > boundaryDicts(boundary.size())
Definition: createBlockMesh.H:64
Foam::faceListList
List< faceList > faceListList
A List of faceList.
Definition: faceListFwd.H:49
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
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:294
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78