Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
39Foam::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
68Foam::labelList Foam::polyMesh::facePatchFaceCells
69(
70 const faceList& patchFaces,
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
126void 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,
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 }
433
434 // Reset the size of the face list
435 faces_.setSize(nFaces);
436}
437
438
439Foam::polyMesh::polyMesh
440(
441 const IOobject& io,
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:
453 objectRegistry(io),
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,
490 IOobject::NO_READ,
491 IOobject::AUTO_WRITE
492 ),
493 Foam::zero{}
494 ),
495 neighbour_
496 (
497 IOobject
498 (
499 "neighbour",
500 instance(),
501 meshSubDir,
502 *this,
503 IOobject::NO_READ,
504 IOobject::AUTO_WRITE
505 ),
506 Foam::zero{}
507 ),
508 clearedPrimitives_(false),
509 boundary_
510 (
511 IOobject
512 (
513 "boundary",
514 instance(),
515 meshSubDir,
516 *this,
517 IOobject::NO_READ,
518 IOobject::AUTO_WRITE
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 pointZones_
528 (
529 IOobject
530 (
531 "pointZones",
532 instance(),
533 meshSubDir,
534 *this,
535 IOobject::NO_READ,
536 IOobject::NO_WRITE
537 ),
538 *this,
539 Foam::zero{}
540 ),
541 faceZones_
542 (
543 IOobject
544 (
545 "faceZones",
546 instance(),
547 meshSubDir,
548 *this,
549 IOobject::NO_READ,
550 IOobject::NO_WRITE
551 ),
552 *this,
553 Foam::zero{}
554 ),
555 cellZones_
556 (
557 IOobject
558 (
559 "cellZones",
560 instance(),
561 meshSubDir,
562 *this,
563 IOobject::NO_READ,
564 IOobject::NO_WRITE
565 ),
566 *this,
567 Foam::zero{}
568 ),
569 moving_(false),
570 topoChanging_(false),
571 storeOldCellCentres_(false),
572 curMotionTimeIndex_(time().timeIndex())
573{
575 << "Constructing polyMesh from cell and boundary shapes." << endl;
576
577 // Calculate faces and cells
578 labelList patchSizes;
579 labelList patchStarts;
580 label defaultPatchStart;
581 label nFaces;
583 setTopology
584 (
585 cellsAsShapes,
586 boundaryFaces,
587 boundaryPatchNames,
588 patchSizes,
589 patchStarts,
590 defaultPatchStart,
591 nFaces,
592 cells
593 );
594
595 // Warning: Patches can only be added once the face list is
596 // completed, as they hold a subList of the face list
597 forAll(boundaryFaces, patchi)
598 {
599 // Add the patch to the list
600 boundary_.set
601 (
602 patchi,
604 (
605 boundaryPatchTypes[patchi],
606 boundaryPatchNames[patchi],
607 patchSizes[patchi],
608 patchStarts[patchi],
609 patchi,
610 boundary_
611 )
612 );
613
614 if
615 (
616 boundaryPatchPhysicalTypes.size()
617 && boundaryPatchPhysicalTypes[patchi].size()
618 )
619 {
620 boundary_[patchi].physicalType() =
621 boundaryPatchPhysicalTypes[patchi];
622 }
623 }
624
625 label nAllPatches = boundaryFaces.size();
626
627 label nDefaultFaces = nFaces - defaultPatchStart;
628 if (syncPar)
629 {
630 reduce(nDefaultFaces, sumOp<label>());
631 }
632
633 if (nDefaultFaces > 0)
634 {
636 << "Found " << nDefaultFaces
637 << " undefined faces in mesh; adding to default patch "
638 << defaultBoundaryPatchName << endl;
639
640 // Check if there already exists a defaultFaces patch as last patch
641 // and reuse it.
642 label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
643
644 if (patchi != -1)
645 {
646 if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
647 {
649 << "Default patch " << boundary_[patchi].name()
650 << " already has faces in it or is not"
651 << " last in list of patches." << exit(FatalError);
652 }
653
655 << "Reusing existing patch " << patchi
656 << " for undefined faces." << endl;
657
658 boundary_.set
659 (
660 patchi,
662 (
663 boundaryPatchTypes[patchi],
664 boundaryPatchNames[patchi],
665 nFaces - defaultPatchStart,
666 defaultPatchStart,
667 patchi,
668 boundary_
669 )
670 );
671 }
672 else
673 {
674 boundary_.set
675 (
676 nAllPatches,
678 (
679 defaultBoundaryPatchType,
680 defaultBoundaryPatchName,
681 nFaces - defaultPatchStart,
682 defaultPatchStart,
683 boundary_.size() - 1,
684 boundary_
685 )
686 );
687
688 nAllPatches++;
689 }
690 }
691
692 // Reset the size of the boundary
693 boundary_.setSize(nAllPatches);
694
695 // Set the primitive mesh
696 initMesh(cells);
697
698 if (syncPar)
699 {
700 // Calculate topology for the patches (processor-processor comms etc.)
701 boundary_.updateMesh();
702
703 // Calculate the geometry for the patches (transformation tensors etc.)
704 boundary_.calcGeometry();
705 }
706
707 if (debug)
708 {
709 if (checkMesh())
711 Info<< "Mesh OK" << endl;
712 }
713 }
714}
715
716
717Foam::polyMesh::polyMesh
718(
719 const IOobject& io,
721 const cellShapeList& cellsAsShapes,
722 const faceListList& boundaryFaces,
723 const wordList& boundaryPatchNames,
725 const word& defaultBoundaryPatchName,
726 const word& defaultBoundaryPatchType,
727 const bool syncPar
728)
729:
732 data_(static_cast<const objectRegistry&>(*this)),
733 points_
734 (
736 (
737 "points",
738 instance(),
739 meshSubDir,
740 *this,
741 IOobject::NO_READ,
742 IOobject::AUTO_WRITE
743 ),
744 std::move(points)
745 ),
746 faces_
747 (
749 (
750 "faces",
751 instance(),
752 meshSubDir,
753 *this,
754 IOobject::NO_READ,
755 IOobject::AUTO_WRITE
756 ),
757 Foam::zero{}
758 ),
759 owner_
760 (
762 (
763 "owner",
764 instance(),
765 meshSubDir,
766 *this,
767 IOobject::NO_READ,
768 IOobject::AUTO_WRITE
769 ),
770 Foam::zero{}
771 ),
772 neighbour_
773 (
775 (
776 "neighbour",
777 instance(),
778 meshSubDir,
779 *this,
780 IOobject::NO_READ,
781 IOobject::AUTO_WRITE
782 ),
783 Foam::zero{}
784 ),
785 clearedPrimitives_(false),
786 boundary_
787 (
789 (
790 "boundary",
791 instance(),
792 meshSubDir,
793 *this,
794 IOobject::NO_READ,
795 IOobject::AUTO_WRITE
796 ),
797 *this,
798 boundaryFaces.size() + 1 // Add room for a default patch
799 ),
800 bounds_(points_, syncPar),
801 comm_(UPstream::worldComm),
802 geometricD_(Zero),
803 solutionD_(Zero),
804 pointZones_
805 (
807 (
808 "pointZones",
809 instance(),
810 meshSubDir,
811 *this,
812 IOobject::NO_READ,
813 IOobject::NO_WRITE
814 ),
815 *this,
816 Foam::zero{}
817 ),
818 faceZones_
819 (
821 (
822 "faceZones",
823 instance(),
824 meshSubDir,
825 *this,
826 IOobject::NO_READ,
827 IOobject::NO_WRITE
828 ),
829 *this,
830 Foam::zero{}
831 ),
832 cellZones_
833 (
835 (
836 "cellZones",
837 instance(),
838 meshSubDir,
839 *this,
840 IOobject::NO_READ,
841 IOobject::NO_WRITE
842 ),
843 *this,
844 Foam::zero{}
845 ),
846 moving_(false),
847 topoChanging_(false),
848 storeOldCellCentres_(false),
849 curMotionTimeIndex_(time().timeIndex())
850{
852 << "Constructing polyMesh from cell and boundary shapes." << endl;
853
854 // Calculate faces and cells
855 labelList patchSizes;
856 labelList patchStarts;
857 label defaultPatchStart;
858 label nFaces;
860 setTopology
861 (
862 cellsAsShapes,
863 boundaryFaces,
864 boundaryPatchNames,
865 patchSizes,
866 patchStarts,
867 defaultPatchStart,
868 nFaces,
869 cells
870 );
871
872 // Warning: Patches can only be added once the face list is
873 // completed, as they hold a subList of the face list
874 forAll(boundaryDicts, patchi)
875 {
876 dictionary patchDict(boundaryDicts[patchi]);
877
878 patchDict.set("nFaces", patchSizes[patchi]);
879 patchDict.set("startFace", patchStarts[patchi]);
880
881 // Add the patch to the list
882 boundary_.set
883 (
884 patchi,
886 (
887 boundaryPatchNames[patchi],
888 patchDict,
889 patchi,
890 boundary_
891 )
892 );
893 }
894
895 label nAllPatches = boundaryFaces.size();
896
897 label nDefaultFaces = nFaces - defaultPatchStart;
898 if (syncPar)
899 {
900 reduce(nDefaultFaces, sumOp<label>());
901 }
902
903 if (nDefaultFaces > 0)
904 {
906 << "Found " << nDefaultFaces
907 << " undefined faces in mesh; adding to default patch "
908 << defaultBoundaryPatchName << endl;
909
910 // Check if there already exists a defaultFaces patch as last patch
911 // and reuse it.
912 label patchi = boundaryPatchNames.find(defaultBoundaryPatchName);
913
914 if (patchi != -1)
915 {
916 if (patchi != boundaryFaces.size()-1 || boundary_[patchi].size())
917 {
919 << "Default patch " << boundary_[patchi].name()
920 << " already has faces in it or is not"
921 << " last in list of patches." << exit(FatalError);
922 }
923
925 << "Reusing existing patch " << patchi
926 << " for undefined faces." << endl;
927
928 boundary_.set
929 (
930 patchi,
932 (
933 boundary_[patchi].type(),
934 boundary_[patchi].name(),
935 nFaces - defaultPatchStart,
936 defaultPatchStart,
937 patchi,
938 boundary_
939 )
940 );
941 }
942 else
943 {
944 boundary_.set
945 (
946 nAllPatches,
948 (
949 defaultBoundaryPatchType,
950 defaultBoundaryPatchName,
951 nFaces - defaultPatchStart,
952 defaultPatchStart,
953 boundary_.size() - 1,
954 boundary_
955 )
956 );
957
958 nAllPatches++;
959 }
960 }
961
962 // Reset the size of the boundary
963 boundary_.setSize(nAllPatches);
964
965 // Set the primitive mesh
966 initMesh(cells);
967
968 if (syncPar)
969 {
970 // Calculate topology for the patches (processor-processor comms etc.)
971 boundary_.updateMesh();
972
973 // Calculate the geometry for the patches (transformation tensors etc.)
974 boundary_.calcGeometry();
975 }
976
977 if (debug)
978 {
979 if (checkMesh())
980 {
981 Info<< "Mesh OK" << endl;
982 }
983 }
984}
985
986
987// ************************************************************************* //
bool found
label size() const noexcept
Definition HashTable.H:358
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Inter-processor communications stream.
Definition UPstream.H:69
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
static bool sameVertices(const face &a, const face &b)
True if the faces have all the same vertices.
Definition face.C:374
Registry of regIOobjects.
const Time & time() const noexcept
Return time registry.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
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.
Cell-face mesh analysis engine.
virtual bool checkMesh(const bool report=false) const
Check mesh for correctness. Returns false for no error.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
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())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
const cellShapeList & cells
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
List< faceList > faceListList
List of faceList.
Definition faceListFwd.H:44
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< cellShape > cellShapeList
List of cellShape.
label timeIndex
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299