Loading...
Searching...
No Matches
polyMesh.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-2017, 2020 OpenFOAM Foundation
9 Copyright (C) 2016-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 "cellIOList.H"
32#include "wedgePolyPatch.H"
33#include "emptyPolyPatch.H"
34#include "globalMeshData.H"
35#include "processorPolyPatch.H"
37#include "indexedOctree.H"
38#include "treeDataCell.H"
39#include "MeshObject.H"
40#include "pointMesh.H"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
47}
48
50
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56void Foam::polyMesh::calcDirections() const
57{
58 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
59 {
60 solutionD_[cmpt] = 1;
61 }
62
63 // Knock out empty and wedge directions. Note:they will be present on all
64 // domains.
65
66 bool hasEmptyPatches = false;
67 bool hasWedgePatches = false;
68
69 vector emptyDirVec = Zero;
70 vector wedgeDirVec = Zero;
71
72 forAll(boundaryMesh(), patchi)
73 {
74 const polyPatch& pp = boundaryMesh()[patchi];
76 {
77 // Force calculation of geometric properties, independent of
78 // size. This avoids parallel synchronisation problems.
79 const vectorField::subField fa(pp.faceAreas());
80
81 if (pp.size())
82 {
83 hasEmptyPatches = true;
84 emptyDirVec += sum(cmptMag(fa));
85 }
86 }
87 else if (isA<wedgePolyPatch>(pp))
88 {
89 const wedgePolyPatch& wpp = refCast<const wedgePolyPatch>(pp);
90
91 // Force calculation of geometric properties, independent of
92 // size. This avoids parallel synchronisation problems.
93 (void)wpp.faceNormals();
94
95 if (pp.size())
96 {
97 hasWedgePatches = true;
98 wedgeDirVec += cmptMag(wpp.centreNormal());
99 }
100 }
101 }
102
103
104 if (returnReduceOr(hasEmptyPatches))
105 {
106 reduce(emptyDirVec, sumOp<vector>());
107
108 emptyDirVec.normalise();
109
110 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
111 {
112 if (emptyDirVec[cmpt] > 1e-6)
113 {
114 solutionD_[cmpt] = -1;
115 }
116 else
117 {
118 solutionD_[cmpt] = 1;
119 }
120 }
121 }
122
123
124 // Knock out wedge directions
125
126 geometricD_ = solutionD_;
127
128 if (returnReduceOr(hasWedgePatches))
129 {
130 reduce(wedgeDirVec, sumOp<vector>());
131
132 wedgeDirVec.normalise();
133
134 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
135 {
136 if (wedgeDirVec[cmpt] > 1e-6)
137 {
138 geometricD_[cmpt] = -1;
139 }
140 else
141 {
142 geometricD_[cmpt] = 1;
143 }
144 }
145 }
146}
147
148
149Foam::autoPtr<Foam::labelIOList> Foam::polyMesh::readTetBasePtIs() const
150{
152 (
153 "tetBasePtIs",
154 instance(),
155 meshSubDir,
156 *this,
160 );
161
162 if (io.typeHeaderOk<labelIOList>(true))
163 {
165 }
167 return nullptr;
168}
169
170
171// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
172
173Foam::polyMesh::polyMesh(const IOobject& io, const bool doInit)
174:
177 data_(static_cast<const objectRegistry&>(*this)),
178 points_
179 (
181 (
182 "points",
183 time().findInstance(meshDir(), "points"),
184 meshSubDir,
185 *this,
186 IOobject::MUST_READ,
187 IOobject::NO_WRITE
188 )
189 ),
190 faces_
191 (
193 (
194 "faces",
195 time().findInstance(meshDir(), "faces"),
196 meshSubDir,
197 *this,
198 IOobject::MUST_READ,
199 IOobject::NO_WRITE
200 )
201 ),
202 owner_
203 (
205 (
206 "owner",
207 faces_.instance(),
208 meshSubDir,
209 *this,
210 IOobject::READ_IF_PRESENT,
211 IOobject::NO_WRITE
212 )
213 ),
214 neighbour_
215 (
217 (
218 "neighbour",
219 faces_.instance(),
220 meshSubDir,
221 *this,
222 IOobject::READ_IF_PRESENT,
223 IOobject::NO_WRITE
224 )
225 ),
226 clearedPrimitives_(false),
227 boundary_
228 (
230 (
231 "boundary",
232 time().findInstance // allow 'newer' boundary file
233 (
234 meshDir(),
235 "boundary",
236 IOobject::MUST_READ,
237 faces_.instance()
238 ),
239 meshSubDir,
240 *this,
241 IOobject::MUST_READ,
242 IOobject::NO_WRITE
243 ),
244 *this
245 ),
246 bounds_(points_),
247 comm_(UPstream::worldComm),
248 geometricD_(Zero),
249 solutionD_(Zero),
250 pointZones_
251 (
253 (
254 "pointZones",
255 faces_.instance(),
256 meshSubDir,
257 *this,
258 IOobject::READ_IF_PRESENT,
259 IOobject::NO_WRITE
260 ),
261 *this,
262 PtrList<entry>()
263 ),
264 faceZones_
265 (
267 (
268 "faceZones",
269 faces_.instance(),
270 meshSubDir,
271 *this,
272 IOobject::READ_IF_PRESENT,
273 IOobject::NO_WRITE
274 ),
275 *this,
276 PtrList<entry>()
277 ),
278 cellZones_
279 (
281 (
282 "cellZones",
283 faces_.instance(),
284 meshSubDir,
285 *this,
286 IOobject::READ_IF_PRESENT,
287 IOobject::NO_WRITE
288 ),
289 *this,
290 PtrList<entry>()
291 ),
292 moving_(false),
293 topoChanging_(false),
294 storeOldCellCentres_(false),
295 curMotionTimeIndex_(time().timeIndex())
296{
297 if (owner_.hasHeaderClass())
298 {
299 initMesh();
300 }
301 else
302 {
303 cellCompactIOList cLst
304 (
305 IOobject
306 (
307 "cells",
308 time().findInstance(meshDir(), "cells"),
309 meshSubDir,
310 *this,
311 IOobject::MUST_READ,
312 IOobject::NO_WRITE
313 )
314 );
315
316 // Set the primitive mesh
317 initMesh(cLst);
318
319 owner_.write();
320 neighbour_.write();
321 }
322
323 if (returnReduceOr(boundary_.empty()))
324 {
325 WarningInFunction
326 << "Missing mesh boundary on one or more domains" << endl;
327
328 // Warn if global empty mesh
329 if (returnReduceAnd(!nPoints()))
330 {
331 WarningInFunction
332 << "No points in mesh" << endl;
333 }
335 {
337 << "No cells in mesh" << endl;
338 }
339 }
340
341 if (doInit)
342 {
343 polyMesh::init(false); // do not init lower levels
344 }
345}
346
347
348bool Foam::polyMesh::init(const bool doInit)
349{
350 if (doInit)
351 {
352 primitiveMesh::init(doInit);
353 }
354
355 // Calculate topology for the patches (processor-processor comms etc.)
356 boundary_.updateMesh();
357
358 // Calculate the geometry for the patches (transformation tensors etc.)
359 boundary_.calcGeometry();
360
361 // Initialise demand-driven data
362 calcDirections();
363
364 return false;
365}
366
367
368Foam::polyMesh::polyMesh
369(
370 const IOobject& io,
372 faceList&& faces,
373 labelList&& owner,
374 labelList&& neighbour,
375 const bool syncPar
376)
377:
380 data_(static_cast<const objectRegistry&>(*this)),
381 points_
382 (
384 (
385 "points",
386 instance(),
387 meshSubDir,
388 *this,
389 IOobject::NO_READ, //io.readOpt(),
390 io.writeOpt()
391 ),
392 std::move(points)
393 ),
394 faces_
395 (
397 (
398 "faces",
399 instance(),
400 meshSubDir,
401 *this,
402 IOobject::NO_READ, //io.readOpt(),
403 io.writeOpt()
404 ),
405 std::move(faces)
406 ),
407 owner_
408 (
410 (
411 "owner",
412 instance(),
413 meshSubDir,
414 *this,
415 IOobject::NO_READ, //io.readOpt(),
416 io.writeOpt()
417 ),
418 std::move(owner)
419 ),
420 neighbour_
421 (
423 (
424 "neighbour",
425 instance(),
426 meshSubDir,
427 *this,
428 IOobject::NO_READ, //io.readOpt(),
429 io.writeOpt()
430 ),
431 std::move(neighbour)
432 ),
433 clearedPrimitives_(false),
434 boundary_
435 (
437 (
438 "boundary",
439 instance(),
440 meshSubDir,
441 *this,
442 IOobject::NO_READ, // ignore since no alternative can be supplied
443 io.writeOpt()
444 ),
445 *this,
447 ),
448 bounds_(points_, syncPar),
449 comm_(UPstream::worldComm),
450 geometricD_(Zero),
451 solutionD_(Zero),
452 pointZones_
453 (
455 (
456 "pointZones",
457 instance(),
458 meshSubDir,
459 *this,
460 IOobject::NO_READ, // ignore since no alternative can be supplied
461 IOobject::NO_WRITE
462 ),
463 *this,
464 Foam::zero{}
465 ),
466 faceZones_
467 (
469 (
470 "faceZones",
471 instance(),
472 meshSubDir,
473 *this,
474 IOobject::NO_READ, // ignore since no alternative can be supplied
475 IOobject::NO_WRITE
476 ),
477 *this,
478 Foam::zero{}
479 ),
480 cellZones_
481 (
483 (
484 "cellZones",
485 instance(),
486 meshSubDir,
487 *this,
488 IOobject::NO_READ, // ignore since no alternative can be supplied
489 IOobject::NO_WRITE
490 ),
491 *this,
492 Foam::zero{}
493 ),
494 moving_(false),
495 topoChanging_(false),
496 storeOldCellCentres_(false),
497 curMotionTimeIndex_(time().timeIndex())
498{
499 // Check if the faces and cells are valid
500 forAll(faces_, facei)
501 {
502 const face& curFace = faces_[facei];
503
504 if (min(curFace) < 0 || max(curFace) > points_.size())
505 {
507 << "Face " << facei << "contains vertex labels out of range: "
508 << curFace << " Max point index = " << points_.size()
509 << abort(FatalError);
510 }
512
513 // Set the primitive mesh
514 initMesh();
515}
516
517
518Foam::polyMesh::polyMesh
519(
520 const IOobject& io,
522 faceList&& faces,
523 cellList&& cells,
524 const bool syncPar
525)
526:
527 objectRegistry(io),
528 primitiveMesh(),
529 data_(static_cast<const objectRegistry&>(*this)),
530 points_
531 (
532 IOobject
533 (
534 "points",
535 instance(),
536 meshSubDir,
537 *this,
538 IOobject::NO_READ,
539 io.writeOpt()
540 ),
541 std::move(points)
542 ),
543 faces_
544 (
545 IOobject
546 (
547 "faces",
548 instance(),
549 meshSubDir,
550 *this,
551 IOobject::NO_READ,
552 io.writeOpt()
553 ),
554 std::move(faces)
555 ),
556 owner_
557 (
558 IOobject
559 (
560 "owner",
561 instance(),
562 meshSubDir,
563 *this,
564 IOobject::NO_READ,
565 io.writeOpt()
566 ),
567 Foam::zero{}
568 ),
569 neighbour_
570 (
571 IOobject
572 (
573 "neighbour",
574 instance(),
575 meshSubDir,
576 *this,
577 IOobject::NO_READ,
578 io.writeOpt()
579 ),
580 Foam::zero{}
581 ),
582 clearedPrimitives_(false),
583 boundary_
584 (
585 IOobject
586 (
587 "boundary",
588 instance(),
589 meshSubDir,
590 *this,
591 IOobject::NO_READ,
592 io.writeOpt()
593 ),
594 *this,
595 Foam::zero{}
596 ),
597 bounds_(points_, syncPar),
598 comm_(UPstream::worldComm),
599 geometricD_(Zero),
600 solutionD_(Zero),
601 pointZones_
602 (
603 IOobject
604 (
605 "pointZones",
606 instance(),
607 meshSubDir,
608 *this,
609 IOobject::NO_READ,
610 IOobject::NO_WRITE
611 ),
612 *this,
613 Foam::zero{}
614 ),
615 faceZones_
616 (
617 IOobject
618 (
619 "faceZones",
620 instance(),
621 meshSubDir,
622 *this,
623 IOobject::NO_READ,
624 IOobject::NO_WRITE
625 ),
626 *this,
627 Foam::zero{}
628 ),
629 cellZones_
630 (
631 IOobject
632 (
633 "cellZones",
634 instance(),
635 meshSubDir,
636 *this,
637 IOobject::NO_READ,
638 IOobject::NO_WRITE
639 ),
640 *this,
641 Foam::zero{}
642 ),
643 moving_(false),
644 topoChanging_(false),
645 storeOldCellCentres_(false),
646 curMotionTimeIndex_(time().timeIndex())
647{
648 // Check if faces are valid
649 forAll(faces_, facei)
650 {
651 const face& curFace = faces_[facei];
652
653 if (min(curFace) < 0 || max(curFace) > points_.size())
654 {
656 << "Face " << facei << "contains vertex labels out of range: "
657 << curFace << " Max point index = " << points_.size()
658 << abort(FatalError);
659 }
660 }
661
662 // Transfer in cell list
663 cellList cLst(std::move(cells));
664
665 // Check if cells are valid
666 forAll(cLst, celli)
667 {
668 const cell& curCell = cLst[celli];
669
670 if (min(curCell) < 0 || max(curCell) > faces_.size())
671 {
673 << "Cell " << celli << "contains face labels out of range: "
674 << curCell << " Max face index = " << faces_.size()
675 << abort(FatalError);
676 }
678
679 // Set the primitive mesh
680 initMesh(cLst);
681}
682
683
684Foam::polyMesh::polyMesh
685(
686 const IOobject& io,
687 const Foam::zero,
688 const bool syncPar
689)
690:
691 polyMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
692{}
693
694
696(
698 autoPtr<faceList>&& faces,
699 autoPtr<labelList>&& owner,
700 autoPtr<labelList>&& neighbour,
701 const labelUList& patchSizes,
702 const labelUList& patchStarts,
703 const bool validBoundary
704)
705{
706 // Clear addressing. Keep geometric props and updateable props for mapping.
707 clearAddressing(true);
708
709 // Take over new primitive data.
710 // Optimized to avoid overwriting data at all
711 if (points)
712 {
713 points_.transfer(*points);
714 bounds_ = boundBox(points_, validBoundary);
715 }
716
717 if (faces)
718 {
719 faces_.transfer(*faces);
720 }
721
722 if (owner)
723 {
724 owner_.transfer(*owner);
725 }
726
727 if (neighbour)
728 {
729 neighbour_.transfer(*neighbour);
730 }
731
732
733 // Reset patch sizes and starts
734 forAll(boundary_, patchi)
735 {
736 boundary_[patchi] = polyPatch
737 (
738 boundary_[patchi],
739 boundary_,
740 patchi,
741 patchSizes[patchi],
742 patchStarts[patchi]
743 );
744 }
745
746
747 // Flags the mesh files as being changed
748 setInstance(time().timeName());
749
750 // Check if the faces and cells are valid
751 forAll(faces_, facei)
752 {
753 const face& curFace = faces_[facei];
754
755 if (min(curFace) < 0 || max(curFace) > points_.size())
756 {
758 << "Face " << facei << " contains vertex labels out of range: "
759 << curFace << " Max point index = " << points_.size()
760 << abort(FatalError);
761 }
762 }
763
764
765 // Set the primitive mesh from the owner_, neighbour_.
766 // Works out from patch end where the active faces stop.
767 initMesh();
768
769
770 if (validBoundary)
771 {
772 // Note that we assume that all the patches stay the same and are
773 // correct etc. so we can already use the patches to do
774 // processor-processor comms.
775
776 // Calculate topology for the patches (processor-processor comms etc.)
777 boundary_.updateMesh();
778
779 // Calculate the geometry for the patches (transformation tensors etc.)
780 boundary_.calcGeometry();
781
782 // Warn if global empty mesh
783 if (returnReduceAnd(!nPoints()) || returnReduceAnd(!nCells()))
784 {
786 << "No points or no cells in mesh" << endl;
788 }
789}
790
791
792// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
793
795{
797 resetMotion();
798}
799
800
801// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
803const Foam::word& Foam::polyMesh::regionName(const word& region)
804{
805 return (region == polyMesh::defaultRegion ? word::null : region);
806}
807
808
810{
811 if (region.empty() || region == polyMesh::defaultRegion)
812 {
814 }
816 return (region / polyMesh::meshSubDir);
817}
818
819
820// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
821
823{
825 {
826 return parent().dbDir();
827 }
828
829 return objectRegistry::dbDir();
830}
831
836}
837
842}
843
846{
847 return points_.instance();
848}
849
852{
853 return faces_.instance();
854}
855
856
858{
859 if (geometricD_.x() == 0)
860 {
861 calcDirections();
862 }
863
864 return geometricD_;
865}
866
868Foam::label Foam::polyMesh::nGeometricD() const
869{
871}
872
873
875{
876 if (solutionD_.x() == 0)
877 {
878 calcDirections();
879 }
880
881 return solutionD_;
882}
883
885Foam::label Foam::polyMesh::nSolutionD() const
886{
888}
889
890
892{
893 if (!tetBasePtIsPtr_)
894 {
895 if (debug)
896 {
898 << "Forcing storage of base points."
899 << endl;
900 }
901
902 labelList basePts
903 (
905 );
906
907 tetBasePtIsPtr_.reset
908 (
909 new labelIOList
910 (
911 IOobject
912 (
913 "tetBasePtIs",
914 instance(),
915 meshSubDir,
916 *this,
920 ),
921 std::move(basePts)
922 )
923 );
925
926 return *tetBasePtIsPtr_;
927}
928
929
932{
933 if (!cellTreePtr_)
934 {
935 Random rndGen(261782);
936
937 treeBoundBox overallBb(points());
938 overallBb.inflate(rndGen, 1e-4, ROOTVSMALL);
939
940 cellTreePtr_.reset
941 (
942 new indexedOctree<treeDataCell>
943 (
944 treeDataCell
945 (
946 false, // not cache bb
947 *this,
948 CELL_TETS // use tet-decomposition for any inside test
949 ),
950 overallBb,
951 8, // maxLevel
952 10, // leafsize
953 5.0 // duplicity
954 )
955 );
956 }
957
958 return *cellTreePtr_;
959}
960
961
963(
964 polyPatchList& plist,
965 const bool validBoundary
966)
967{
968 if (boundaryMesh().size())
969 {
971 << "boundary already exists"
972 << abort(FatalError);
973 }
974
975 // Reset valid directions
976 geometricD_ = Zero;
977 solutionD_ = Zero;
978
979 boundary_.transfer(plist);
980
981 // parallelData depends on the processorPatch ordering so force
982 // recalculation. Problem: should really be done in removeBoundary but
983 // there is some info in parallelData which might be interesting inbetween
984 // removeBoundary and addPatches.
985 globalMeshDataPtr_.reset(nullptr);
986
987 if (validBoundary)
988 {
989 // Calculate topology for the patches (processor-processor comms etc.)
990 boundary_.updateMesh();
991
992 // Calculate the geometry for the patches (transformation tensors etc.)
993 boundary_.calcGeometry();
994
995 boundary_.checkDefinition();
996 }
997}
998
999
1001(
1002 PtrList<pointZone>&& pz,
1003 PtrList<faceZone>&& fz,
1004 PtrList<cellZone>&& cz
1005)
1006{
1007 if (pointZones_.size() || faceZones_.size() || cellZones_.size())
1008 {
1010 << "point, face or cell zone already exists"
1011 << abort(FatalError);
1012 }
1013
1014 // Point zones - take ownership of the pointers
1015 if (pz.size())
1016 {
1017 pointZones_.clear();
1018 pointZones_.transfer(pz);
1019 pointZones_.writeOpt(IOobject::AUTO_WRITE);
1020 }
1021
1022 // Face zones - take ownership of the pointers
1023 if (fz.size())
1024 {
1025 faceZones_.clear();
1026 faceZones_.transfer(fz);
1027 faceZones_.writeOpt(IOobject::AUTO_WRITE);
1028 }
1029
1030 // Cell zones - take ownership of the pointers
1031 if (cz.size())
1032 {
1033 cellZones_.clear();
1034 cellZones_.transfer(cz);
1035 cellZones_.writeOpt(IOobject::AUTO_WRITE);
1036 }
1037}
1038
1039
1041(
1042 const List<polyPatch*>& p,
1043 const bool validBoundary
1044)
1045{
1046 // Acquire ownership of the pointers
1047 polyPatchList plist(const_cast<List<polyPatch*>&>(p));
1048
1049 addPatches(plist, validBoundary);
1050}
1051
1052
1054(
1055 const List<pointZone*>& pz,
1056 const List<faceZone*>& fz,
1057 const List<cellZone*>& cz
1058)
1059{
1060 // Acquire ownership of the pointers
1061 addZones
1062 (
1064 PtrList<faceZone>(const_cast<List<faceZone*>&>(fz)),
1065 PtrList<cellZone>(const_cast<List<cellZone*>&>(cz))
1066 );
1067}
1068
1069
1071{
1072 if (clearedPrimitives_)
1073 {
1075 << "points deallocated"
1077 }
1078
1079 return points_;
1080}
1081
1084{
1085 return io.upToDate(points_);
1086}
1087
1090{
1091 io.eventNo() = points_.eventNo()+1;
1092}
1093
1094
1096{
1097 if (clearedPrimitives_)
1098 {
1100 << "faces deallocated"
1102 }
1103
1104 return faces_;
1105}
1106
1109{
1110 return owner_;
1111}
1112
1115{
1116 return neighbour_;
1117}
1118
1119
1121{
1122 if (!moving_)
1123 {
1124 return points_;
1125 }
1126
1127 if (!oldPointsPtr_)
1128 {
1129 if (debug)
1130 {
1132 }
1133
1134 oldPointsPtr_.reset(new pointField(points_));
1135 curMotionTimeIndex_ = time().timeIndex();
1136 }
1137
1138 return *oldPointsPtr_;
1139}
1140
1141
1143{
1144 storeOldCellCentres_ = true;
1145
1146 if (!moving_)
1147 {
1148 return cellCentres();
1149 }
1150
1151 if (!oldCellCentresPtr_)
1152 {
1153 oldCellCentresPtr_.reset(new pointField(cellCentres()));
1154 }
1155
1156 return *oldCellCentresPtr_;
1157}
1158
1159
1160void Foam::polyMesh::movePoints(const pointField& newPoints)
1161{
1163 << "Moving points for time " << time().value()
1164 << " index " << time().timeIndex() << endl;
1165
1166 if (newPoints.size() != points_.size())
1167 {
1169 << "Size of newPoints " << newPoints.size()
1170 << " does not correspond to current mesh points size "
1171 << points_.size()
1172 << exit(FatalError);
1173 }
1174
1175
1176 moving(true);
1177
1178 // Pick up old points
1179 if (curMotionTimeIndex_ != time().timeIndex())
1180 {
1181 if (debug)
1182 {
1183 Info<< "void polyMesh::movePoints(const pointField&) : "
1184 << " Storing current points for time " << time().value()
1185 << " index " << time().timeIndex() << endl;
1186 }
1187
1188 if (storeOldCellCentres_)
1189 {
1190 oldCellCentresPtr_.reset(nullptr);
1191 oldCellCentresPtr_.reset(new pointField(cellCentres()));
1192 }
1193
1194 // Mesh motion in the new time step
1195 oldPointsPtr_.reset(nullptr);
1196 oldPointsPtr_.reset(new pointField(points_));
1197 curMotionTimeIndex_ = time().timeIndex();
1198 }
1199
1200 points_ = newPoints;
1201
1202 bool moveError = false;
1203 if (debug)
1204 {
1205 // Check mesh motion
1206 if (checkMeshMotion(points_, true))
1207 {
1208 moveError = true;
1209
1211 << "Moving the mesh with given points will "
1212 << "invalidate the mesh." << nl
1213 << "Mesh motion should not be executed." << endl;
1214 }
1215 }
1216
1217 points_.writeOpt(IOobject::AUTO_WRITE);
1218 points_.instance() = time().timeName();
1219 points_.eventNo() = getEvent();
1220
1221 if (tetBasePtIsPtr_)
1222 {
1223 tetBasePtIsPtr_->writeOpt(IOobject::AUTO_WRITE);
1224 tetBasePtIsPtr_->instance() = time().timeName();
1225 tetBasePtIsPtr_->eventNo() = getEvent();
1226 }
1227
1228 // Currently a no-op; earlier versions set meshPhi and call
1229 // primitiveMesh::clearGeom
1230 (void)primitiveMesh::movePoints(points_, oldPoints());
1231
1232 // Update the mesh geometry (via fvGeometryScheme)
1233 // - updateGeom is virtual -> calls fvMesh::updateGeom (or higher)
1234 // - fvMesh::updateGeom defers to surfaceInterpolation::updateGeom(),
1235 // which defers to fvGeometryScheme::movePoints()
1236 // - set the mesh flux
1237 // - clear out/recalculate stale geometry
1238 updateGeom();
1239
1240 // Adjust parallel shared points
1241 if (globalMeshDataPtr_)
1242 {
1243 globalMeshDataPtr_->movePoints(points_);
1244 }
1245
1246 // Force recalculation of all geometric data with new points
1247
1248 bounds_ = boundBox(points_);
1249 boundary_.movePoints(points_);
1250
1251 pointZones_.movePoints(points_);
1252 faceZones_.movePoints(points_);
1253 cellZones_.movePoints(points_);
1254
1255 // Reset cell tree - it gets built from mesh geometry so might have
1256 // wrong boxes. It is correct as long as none of the cells leaves
1257 // the boxes it is in which most likely is almost never the case except
1258 // for tiny displacements. An alternative is to check the displacements
1259 // to see if they are tiny - imagine a big windtunnel with a small rotating
1260 // object. In this case the processors without the rotating object wouldn't
1261 // have to clear any geometry. However your critical path still stays the
1262 // same so no time would be gained (unless the decomposition gets weighted).
1263 // Small benefit for lots of scope for problems so not done.
1264 cellTreePtr_.reset(nullptr);
1265
1266 // Reset valid directions (could change with rotation)
1267 geometricD_ = Zero;
1268 solutionD_ = Zero;
1269
1270 // Note: tet-base decomposition does not get cleared. Ideally your face
1271 // decomposition should not change during mesh motion ...
1272
1273
1276
1277 const_cast<Time&>(time()).functionObjects().movePoints(*this);
1278
1279
1280 if (debug && moveError)
1281 {
1282 // Write mesh to ease debugging. Note we want to avoid calling
1283 // e.g. fvMesh::write since meshPhi not yet complete.
1285 }
1286}
1287
1288
1289void Foam::polyMesh::resetMotion() const
1291 curMotionTimeIndex_ = 0;
1292 oldPointsPtr_.reset(nullptr);
1293 oldCellCentresPtr_.reset(nullptr);
1294}
1295
1298{
1299 return bool(globalMeshDataPtr_);
1300}
1301
1302
1304{
1305 if (!globalMeshDataPtr_)
1306 {
1307 if (debug)
1308 {
1309 Pout<< "polyMesh::globalData() const : "
1310 << "Constructing parallelData from processor topology"
1311 << endl;
1312 }
1313 // Construct globalMeshData using processorPatch information only.
1314 globalMeshDataPtr_.reset(new globalMeshData(*this));
1315 }
1316
1317 return *globalMeshDataPtr_;
1318}
1319
1320
1321void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1322{
1323 fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1324
1325 rm(meshFilesPath/"points");
1326 rm(meshFilesPath/"faces");
1327 rm(meshFilesPath/"owner");
1328 rm(meshFilesPath/"neighbour");
1329 rm(meshFilesPath/"cells");
1330 rm(meshFilesPath/"boundary");
1331 rm(meshFilesPath/"pointZones");
1332 rm(meshFilesPath/"faceZones");
1333 rm(meshFilesPath/"cellZones");
1334 rm(meshFilesPath/"meshModifiers");
1335 rm(meshFilesPath/"parallelData");
1336
1337 // remove subdirectories
1338 if (isDir(meshFilesPath/"sets"))
1339 {
1340 rmDir(meshFilesPath/"sets");
1341 }
1342}
1343
1345void Foam::polyMesh::removeFiles() const
1346{
1348}
1349
1350
1352(
1353 const point& p,
1354 label& celli,
1355 label& tetFacei,
1356 label& tetPti
1357) const
1358{
1359 celli = -1;
1360 tetFacei = -1;
1361 tetPti = -1;
1362
1363 const indexedOctree<treeDataCell>& tree = cellTree();
1364
1365 // Find point inside cell
1366 celli = tree.findInside(p);
1367
1368 if (celli != -1)
1370 // Check the nearest cell to see if the point is inside.
1371 findTetFacePt(celli, p, tetFacei, tetPti);
1372 }
1373}
1374
1375
1377(
1378 const label celli,
1379 const point& p,
1380 label& tetFacei,
1381 label& tetPti
1382) const
1383{
1384 const polyMesh& mesh = *this;
1387 tetFacei = tet.face();
1388 tetPti = tet.tetPt();
1389}
1390
1391
1393(
1394 const point& p,
1395 label celli,
1396 const cellDecomposition decompMode
1397) const
1398{
1399 switch (decompMode)
1400 {
1401 case FACE_PLANES:
1402 {
1403 return primitiveMesh::pointInCell(p, celli);
1404 }
1405 break;
1406
1407 case FACE_CENTRE_TRIS:
1408 {
1409 // only test that point is on inside of plane defined by cell face
1410 // triangles
1411 const cell& cFaces = cells()[celli];
1412
1413 forAll(cFaces, cFacei)
1414 {
1415 label facei = cFaces[cFacei];
1416 const face& f = faces_[facei];
1417 const point& fc = faceCentres()[facei];
1418 bool isOwn = (owner_[facei] == celli);
1419
1420 forAll(f, fp)
1421 {
1422 label pointi;
1423 label nextPointi;
1424
1425 if (isOwn)
1426 {
1427 pointi = f[fp];
1428 nextPointi = f.nextLabel(fp);
1429 }
1430 else
1431 {
1432 pointi = f.nextLabel(fp);
1433 nextPointi = f[fp];
1434 }
1435
1436 triPointRef faceTri
1437 (
1438 points()[pointi],
1439 points()[nextPointi],
1440 fc
1441 );
1442
1443 vector proj = p - faceTri.centre();
1444
1445 if ((faceTri.areaNormal() & proj) > 0)
1446 {
1447 return false;
1448 }
1449 }
1450 }
1451 return true;
1452 }
1453 break;
1454
1455 case FACE_DIAG_TRIS:
1456 {
1457 // only test that point is on inside of plane defined by cell face
1458 // triangles
1459 const cell& cFaces = cells()[celli];
1460
1461 forAll(cFaces, cFacei)
1462 {
1463 label facei = cFaces[cFacei];
1464 const face& f = faces_[facei];
1465
1466 for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1467 {
1468 // Get tetIndices of face triangle
1469 tetIndices faceTetIs(celli, facei, tetPti);
1470
1471 triPointRef faceTri = faceTetIs.faceTri(*this);
1472
1473 vector proj = p - faceTri.centre();
1474
1475 if ((faceTri.areaNormal() & proj) > 0)
1476 {
1477 return false;
1478 }
1479 }
1480 }
1481
1482 return true;
1483 }
1484 break;
1485
1486 case CELL_TETS:
1487 {
1488 label tetFacei;
1489 label tetPti;
1490
1491 findTetFacePt(celli, p, tetFacei, tetPti);
1492
1493 return tetFacei != -1;
1494 }
1495 break;
1496 }
1497
1498 return false;
1499}
1500
1501
1502Foam::label Foam::polyMesh::findCell
1503(
1504 const point& p,
1505 const cellDecomposition decompMode
1506) const
1507{
1508 if
1509 (
1511 && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1512 )
1513 {
1514 // Force construction of face-diagonal decomposition before testing
1515 // for zero cells.
1516 //
1517 // If parallel running a local domain might have zero cells so never
1518 // construct the face-diagonal decomposition which uses parallel
1519 // transfers.
1520 (void)tetBasePtIs();
1521 }
1522
1523 if (nCells() == 0)
1524 {
1525 return -1;
1526 }
1527
1528 if (decompMode == CELL_TETS)
1529 {
1530 // Advanced search method utilizing an octree
1531 // and tet-decomposition of the cells
1532
1533 label celli;
1534 label tetFacei;
1535 label tetPti;
1536
1537 findCellFacePt(p, celli, tetFacei, tetPti);
1538
1539 return celli;
1540 }
1541 else
1542 {
1543 // Approximate search avoiding the construction of an octree
1544 // and cell decomposition
1545
1546 if (Pstream::parRun() && decompMode == FACE_DIAG_TRIS)
1547 {
1548 // Force construction of face-diagonal decomposition before testing
1549 // for zero cells. If parallel running a local domain might have
1550 // zero cells so never construct the face-diagonal decomposition
1551 // (which uses parallel transfers)
1552 (void)tetBasePtIs();
1553 }
1554
1555 // Find the nearest cell centre to this location
1556 label celli = findNearestCell(p);
1557
1558 // If point is in the nearest cell return
1559 if (pointInCell(p, celli, decompMode))
1560 {
1561 return celli;
1562 }
1563 else
1564 {
1565 // Point is not in the nearest cell so search all cells
1566
1567 for (label celli = 0; celli < nCells(); celli++)
1568 {
1569 if (pointInCell(p, celli, decompMode))
1570 {
1571 return celli;
1572 }
1573 }
1574
1575 return -1;
1577 }
1578}
1579
1580
1581// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1582
1584(
1585 IOstreamOption streamOpt,
1586 const bool writeOnProc
1587) const
1588{
1589 // Currently no special treatment. Just write the objects
1590
1591 return objectRegistry::writeObject(streamOpt, writeOnProc);
1592}
1593
1594
1595// ************************************************************************* //
if(patchID !=-1)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
SubField< vector > subField
Definition Field.H:183
label size() const noexcept
Definition HashTable.H:358
@ NO_REGISTER
Do not request registration (bool: false).
writeOption writeOpt() const noexcept
Get the write option.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
@ 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
bool hasHeaderClass() const noexcept
True if headerClassName() is non-empty (after reading).
Definition IOobjectI.H:261
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
A simple container for options an IOstream can normally have.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Random number generator.
Definition Random.H:56
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Inter-processor communications stream.
Definition UPstream.H:69
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
bool empty() const noexcept
True if the list is empty (ie, size() is zero).
Definition UPtrListI.H:99
static constexpr direction nComponents
Number of components in this vector space.
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Non-pointer based hierarchical recursive searching.
static void movePoints(objectRegistry &obr)
Update for mesh motion.
Registry of regIOobjects.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the objects using stream options.
virtual const fileName & dbDir() const
Local directory path of this objectRegistry relative to the time.
const objectRegistry & parent() const noexcept
Return the parent objectRegistry.
const Time & time() const noexcept
Return time registry.
label getEvent() const
Return new event number.
static labelList findFaceBasePts(const polyMesh &mesh, scalar tol=minTetQuality, bool report=false)
Find a suitable base point for each face for decomposition into tets.
static tetIndices findTet(const polyMesh &mesh, label cI, const point &pt)
Find the tet decomposition of the cell containing the given point.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual ~polyMesh()
Destructor.
Definition polyMesh.C:787
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition polyMesh.C:1345
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition polyMeshIO.C:29
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition polyMesh.C:861
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition polyMesh.C:1290
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition polyMesh.C:956
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write items held in the objectRegistry. Normally includes mesh components (points,...
Definition polyMesh.C:1577
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir).
Definition polyMesh.C:826
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition polyMesh.H:105
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
bool moving() const noexcept
Is mesh moving.
Definition polyMesh.H:732
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition polyMesh.C:341
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition polyMesh.C:815
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition polyMesh.C:1496
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition polyMesh.C:884
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition polyMesh.C:1370
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition polyMesh.C:1386
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition polyMesh.C:878
virtual const pointField & oldPoints() const
Return old points (mesh motion).
Definition polyMesh.C:1113
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition polyMesh.C:838
void clearAddressing(const bool isMeshUpdate=false)
Clear addressing.
void resetMotion() const
Reset motion.
Definition polyMesh.C:1282
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition polyMesh.C:994
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition polyMesh.C:1076
void removeFiles(const fileName &instanceDir) const
Remove all files from mesh instance.
Definition polyMesh.C:1314
virtual void movePoints(const pointField &)
Move points.
Definition polyMesh.C:1153
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition polyMesh.H:690
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition polyMesh.C:1082
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition polyMesh.C:924
void removeFiles() const
Remove all files from mesh instance().
Definition polyMesh.C:1338
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion).
Definition polyMesh.C:1135
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition polyMesh.C:689
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition polyMesh.C:832
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition polyMesh.C:867
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition polyMesh.C:850
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Cell-face mesh analysis engine.
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
bool pointInCell(const point &p, label celli) const
Return true if the point is in the cell.
const vectorField & faceCentres() const
const vectorField & cellCentres() const
label nCells() const noexcept
Number of mesh cells.
void movePoints(const pointField &p, const pointField &oldP)
Move points.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
const cellList & cells() const
virtual void updateGeom()
Update all geometric data.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
label face() const noexcept
Return the face index.
Definition tetIndices.H:156
triPointRef faceTri(const polyMesh &mesh) const
The triangle geometry for the face for this tet. The normal of the tri points out of the cell.
label tetPt() const noexcept
Return the characterising tet point index.
Definition tetIndices.H:166
Standard boundBox with extra functionality for use in octree.
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
Definition triangleI.H:169
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition triangleI.H:144
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
label nPoints
const cellShapeList & cells
word timeName
Definition getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition POSIX.C:1406
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
bool rmDir(const fileName &directory, const bool silent=false, const bool emptyOnly=false)
Remove a directory and its contents recursively,.
Definition POSIX.C:1435
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const direction noexcept
Definition scalarImpl.H:265
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition POSIX.C:862
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
labelList f(nPoints)
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen