Loading...
Searching...
No Matches
meshRefinement.H
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 OpenFOAM Foundation
9 Copyright (C) 2015-2024 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
27Class
28 Foam::meshRefinement
29
30Description
31 Helper class which maintains intersections of (changing) mesh with
32 (static) surfaces.
33
34 Maintains
35 - per face any intersections of the cc-cc segment with any of the surfaces
36
37SourceFiles
38 meshRefinement.C
39 meshRefinementBaffles.C
40 meshRefinementGapRefine.C
41 meshRefinementMerge.C
42 meshRefinementProblemCells.C
43 meshRefinementRefine.C
44
45\*---------------------------------------------------------------------------*/
46
47#ifndef Foam_meshRefinement_H
48#define Foam_meshRefinement_H
49
50#include "hexRef8.H"
51#include "mapPolyMesh.H"
52#include "autoPtr.H"
53#include "labelPairHashes.H"
55#include "pointFieldsFwd.H"
56#include "Tuple2.H"
57#include "pointIndexHit.H"
58#include "wordPairHashes.H"
59#include "surfaceZonesInfo.H"
60#include "volumeType.H"
61#include "DynamicField.H"
62#include "coordSetWriter.H"
63#include "surfaceWriter.H"
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67namespace Foam
68{
69
70// Forward Declarations
71class fvMesh;
73class mapDistribute;
77class shellSurfaces;
78class removeCells;
80class removePoints;
82class snapParameters;
84/*---------------------------------------------------------------------------*\
85 Class meshRefinement Declaration
86\*---------------------------------------------------------------------------*/
87
88class meshRefinement
89{
90public:
91
92 // Public Data Types
93
94 //- Enumeration for how to operate
96 {
97 CASTELLATED, // split-hex refinement, snapping, layer addition
98 CASTELLATEDBUFFERLAYER, // same but add buffer layers before
99 // snapping
100 CASTELLATEDBUFFERLAYER2 // experimental. wip.
101 };
102
103 static const Enum<MeshType> MeshTypeNames;
104
105
106 //- Enumeration for what to debug. Used as a bit-pattern.
107 enum debugType
109 MESH = (1 << 0),
111 FEATURESEEDS = (1 << 2),
112 ATTRACTION = (1 << 3),
113 LAYERINFO = (1 << 4)
114 };
116 static const Enum<debugType> debugTypeNames;
117
119 //enum outputType
120 //{
121 // OUTPUTLAYERINFO = (1 << 0)
122 //};
123 //
124 //static const Enum<outputType> outputTypeNames;
125
126 //- Enumeration for what to write. Used as a bit-pattern.
127 enum writeType
129 WRITEMESH = (1 << 0),
131 WRITELEVELS = (1 << 2),
132 WRITELAYERSETS = (1 << 3),
134 };
135
136 static const Enum<writeType> writeTypeNames;
138 //- Enumeration for how userdata is to be mapped upon refinement.
139 enum mapType
140 {
141 MASTERONLY = 1,
143 REMOVE = 4
144 };
146 //- Enumeration for what to do with co-planar patch faces on a single
147 // cell
148 enum FaceMergeType
149 {
150 NONE, // no merging
151 GEOMETRIC, // use feature angle
152 IGNOREPATCH // use feature angle, allow merging of different
153 // patches
154 };
155
157private:
159 // Static Data Members
160
161 //- Control of writing level
162 static writeType writeLevel_;
163
165 //static outputType outputLevel_;
166
167
168 // Private Data
169
170 //- Reference to mesh
171 fvMesh& mesh_;
172
173 //- Tolerance used for sorting coordinates (used in 'less' routine)
174 const scalar mergeDistance_;
175
176 //- Overwrite the mesh?
177 const bool overwrite_;
178
179 //- Instance of mesh upon construction. Used when in overwrite_ mode.
180 const word oldInstance_;
181
182 //- All surface-intersection interaction
183 const refinementSurfaces& surfaces_;
184
185 //- All feature-edge interaction
186 const refinementFeatures& features_;
187
188 //- All shell-refinement interaction
189 const shellSurfaces& shells_;
190
191 //- All limit-refinement interaction
192 const shellSurfaces& limitShells_;
193
194 //- How to generate mesh
195 const MeshType meshType_;
196
197 //- Are we operating in test mode?
198 const bool dryRun_;
199
200
201 //- Refinement engine
202 hexRef8 meshCutter_;
203
204 //- Per cc-cc vector the index of the surface hit
205 labelIOList surfaceIndex_;
206
207
208 // For baffle merging
209
210 //- Original patch for baffle faces that used to be on
211 // coupled patches
212 Map<label> faceToCoupledPatch_;
213
214
215 //- User supplied face based data.
216 List<Tuple2<mapType, labelList>> userFaceData_;
217
218 //- Meshed patches - are treated differently. Stored as wordList since
219 // order changes.
220 wordList meshedPatches_;
221
222 //- FaceZone to master patch name
223 HashTable<word> faceZoneToMasterPatch_;
224
225 //- FaceZone to slave patch name
226 HashTable<word> faceZoneToSlavePatch_;
227
228 //- FaceZone to method to handle faces
230
231
232 // Private Member Functions
233
234 //- Add patchfield of given type to all fields on mesh
235 template<class GeoField>
236 static void addPatchFields(fvMesh&, const word& patchFieldType);
237
238 //- Reorder patchfields of all fields on mesh
239 template<class GeoField>
240 static void reorderPatchFields(fvMesh&, const labelList& oldToNew);
241
242 //- Find out which faces have changed given cells (old mesh labels)
243 // that were marked for refinement.
244 static labelList getChangedFaces
245 (
246 const mapPolyMesh&,
247 const labelList& oldCellsToRefine
248 );
249
250 // Debug: count number of master-faces (and do some checking
251 // for consistency)
252 label globalFaceCount(const labelList& elems) const;
253
254 //- Calculate coupled boundary end vector and refinement level
255 void calcNeighbourData(labelList& neiLevel, pointField& neiCc) const;
256
257 //- Calculate rays from cell-centre to cell-centre and corresponding
258 // min cell refinement level
259 void calcCellCellRays
260 (
261 const pointField& neiCc,
262 const labelUList& neiLevel,
263 const labelUList& testFaces,
264 pointField& start,
265 pointField& end,
266 labelList& minLevel
267 ) const;
268
269 //- Get cells which are inside any closed surface. Note that
270 // all closed surfaces
271 // will have already been oriented to have keepPoint outside.
272 labelList getInsideCells(const word&) const;
273
274 //- Do all to remove inside cells
275 autoPtr<mapPolyMesh> removeInsideCells
276 (
277 const string& msg,
278 const label exposedPatchi
279 );
280
281
282 // Refinement candidate selection
283
284 //- Mark cell for refinement (if not already marked). Return false
285 // if refinelimit hit. Keeps running count (in nRefine) of cells
286 // marked for refinement
287 static bool markForRefine
288 (
289 const label markValue,
290 const label nAllowRefine,
291 label& cellValue,
292 label& nRefine
293 );
294
295 //- Mark every cell with level of feature passing through it
296 // (or -1 if not passed through). Uses tracking.
297 void markFeatureCellLevel
298 (
299 const pointField& keepPoints,
300 labelList& maxFeatureLevel
301 ) const;
302
303 //- Calculate list of cells to refine based on intersection of
304 // features.
305 label markFeatureRefinement
306 (
307 const pointField& keepPoints,
308 const label nAllowRefine,
309
311 label& nRefine
312 ) const;
313
314 //- Mark cells for distance-to-feature based refinement.
315 label markInternalDistanceToFeatureRefinement
316 (
317 const label nAllowRefine,
319 label& nRefine
320 ) const;
321
322 //- Mark cells for refinement-shells based refinement.
323 label markInternalRefinement
324 (
325 const label nAllowRefine,
327 label& nRefine
328 ) const;
329
330 //- Unmark cells for refinement based on limit-shells. Return number
331 // of limited cells.
332 label unmarkInternalRefinement
333 (
335 label& nRefine
336 ) const;
337
338 //- Collect faces that are intersected and whose neighbours aren't
339 // yet marked for refinement.
340 labelList getRefineCandidateFaces
341 (
342 const labelList& refineCell
343 ) const;
344
345 //- Mark cells for surface intersection based refinement.
346 label markSurfaceRefinement
347 (
348 const label nAllowRefine,
349 const labelList& neiLevel,
350 const pointField& neiCc,
352 label& nRefine
353 ) const;
354
355 //- Collect cells intersected by the surface that are candidates
356 // for gap checking. Used inside markSurfaceGapRefinement
357 void collectGapCandidates
358 (
359 const shellSurfaces& shells,
360 const labelList& testFaces,
361 const labelList& neiLevel,
362 const pointField& neiCc,
363 labelList& cellToCompact,
364 labelList& bFaceToCompact,
365 List<FixedList<label, 3>>& shellGapInfo,
366 List<volumeType>& shellGapMode
367 ) const;
368 void collectGapCells
369 (
370 const scalar planarCos,
371
372 const List<FixedList<label, 3>>& extendedGapLevel,
373 const List<volumeType>& extendedGapMode,
374 const labelList& testFaces,
375 const pointField& start,
376 const pointField& end,
377
378 const labelList& cellToCompact,
379 const labelList& bFaceToCompact,
380 const List<FixedList<label, 3>>& shellGapInfo,
381 const List<volumeType>& shellGapMode,
382
383 const label nAllowRefine,
384 const labelList& neiLevel,
385 const pointField& neiCc,
386
388 label& nRefine
389 ) const;
390
391
392 //- Mark cells intersected by the surface if they are inside
393 // close gaps
394 label markSurfaceGapRefinement
395 (
396 const scalar planarCos,
397 const label nAllowRefine,
398 const labelList& neiLevel,
399 const pointField& neiCc,
400
402 label& nRefine
403 ) const;
404
405 //- Generate single ray from nearPoint in direction of nearNormal
406 label generateRays
407 (
408 const point& nearPoint,
409 const vector& nearNormal,
410 const FixedList<label, 3>& gapInfo,
411 const volumeType& mode,
412
413 const label cLevel,
414
415 DynamicField<point>& start,
417 ) const;
418
419 //- Generate pairs of rays through cell centre
420 // Each ray pair has start, end, and expected gap size
421 label generateRays
422 (
423 const bool useSurfaceNormal,
424
425 const point& nearPoint,
426 const vector& nearNormal,
427 const FixedList<label, 3>& gapInfo,
428 const volumeType& mode,
429
430 const point& cc,
431 const label cLevel,
432
433 DynamicField<point>& start,
435 DynamicField<scalar>& gapSize,
436
437 DynamicField<point>& start2,
439 DynamicField<scalar>& gapSize2
440 ) const;
441
442 //- Select candidate cells (cells inside a shell with gapLevel
443 // specified)
444 void selectGapCandidates
445 (
446 const labelList& refineCell,
447 const label nRefine,
448
449 labelList& cellMap,
450 labelList& gapShell,
451 List<FixedList<label, 3>>& shellGapInfo,
452 List<volumeType>& shellGapMode
453 ) const;
454
455 //- Merge gap information coming from shell and from surface
456 // (surface wins)
457 void mergeGapInfo
458 (
459 const FixedList<label, 3>& shellGapInfo,
460 const volumeType shellGapMode,
461 const FixedList<label, 3>& surfGapInfo,
462 const volumeType surfGapMode,
463 FixedList<label, 3>& gapInfo,
464 volumeType& gapMode
465 ) const;
466
467 //- Mark cells for non-surface intersection based gap refinement
468 label markInternalGapRefinement
469 (
470 const scalar planarCos,
471 const bool spreadGapSize,
472 const label nAllowRefine,
474 label& nRefine,
475 labelList& numGapCells,
476 scalarField& gapSize
477 ) const;
478
479 //- Refine cells containing small gaps
480 label markSmallFeatureRefinement
481 (
482 const scalar planarCos,
483 const label nAllowRefine,
484 const labelList& neiLevel,
485 const pointField& neiCc,
486
488 label& nRefine
489 ) const;
490
491 //- Refine cells containing triangles with high refinement level
492 // (currently due to high curvature or being inside shell with
493 // high level)
494 label markSurfaceFieldRefinement
495 (
496 const label nAllowRefine,
497 const labelList& neiLevel,
498 const pointField& neiCc,
499
501 label& nRefine
502 ) const;
503
504 //- Helper: count number of normals1 that are in normals2
505 label countMatches
506 (
507 const List<point>& normals1,
508 const List<point>& normals2,
509 const scalar tol = 1e-6
510 ) const;
511
512 //- Detect if two intersection points are high curvature (w.r.t.
513 // lengthscale
514 bool highCurvature
515 (
516 const scalar minCosAngle,
517 const scalar lengthScale,
518 const point& p0,
519 const vector& n0,
520 const point& p1,
521 const vector& n1
522 ) const;
523
524 //- Mark cells for surface curvature based refinement. Marks if
525 // local curvature > curvature or if on different regions
526 // (markDifferingRegions)
527 label markSurfaceCurvatureRefinement
528 (
529 const scalar curvature,
530 const label nAllowRefine,
531 const labelList& neiLevel,
532 const pointField& neiCc,
534 label& nRefine
535 ) const;
536
537 //- Mark cell if local a gap topology or
538 bool checkProximity
539 (
540 const scalar planarCos,
541 const label nAllowRefine,
542
543 const label surfaceLevel,
544 const vector& surfaceLocation,
545 const vector& surfaceNormal,
546
547 const label celli,
548
549 label& cellMaxLevel,
550 vector& cellMaxLocation,
551 vector& cellMaxNormal,
552
554 label& nRefine
555 ) const;
556
557 //- Mark cells for surface proximity based refinement.
558 label markProximityRefinement
559 (
560 const scalar curvature,
561
562 // Per region the min and max cell level
563 const labelList& surfaceMinLevel,
564 const labelList& surfaceMaxLevel,
565
566 const label nAllowRefine,
567 const labelList& neiLevel,
568 const pointField& neiCc,
569
571 label& nRefine
572 ) const;
573
574 //- Mark faces for additional proximity based testing. Extends
575 // testFaces
576 void markProximityCandidateFaces
577 (
578 const labelList& blockedSurfaces,
579 const List<scalarList>& regionToBlockSize,
580 const labelList& neiLevel,
581 const pointField& neiCc,
582
583 labelList& testFaces
584 ) const;
585
586 //- Mark cells for surface proximity based refinement. Uses
587 // ray-shooting from markInternalGapRefinement.
588 label markFakeGapRefinement
589 (
590 const scalar planarCos,
591 const labelList& blockedSurfaces,
592 const label nAllowRefine,
593 const labelList& neiLevel,
594 const pointField& neiCc,
595
597 label& nRefine
598 ) const;
599
600
601 // Baffle handling
602
603 //- Get faces to repatch. Returns map from face to patch.
604 Map<labelPair> getZoneBafflePatches
605 (
606 const bool allowBoundary,
607 const labelList& globalToMasterPatch,
608 const labelList& globalToSlavePatch
609 ) const;
610
611 //- Calculate intersections. Return per face -1 or the global
612 // surface region
613 void getIntersections
614 (
615 const labelList& surfacesToTest,
616 const pointField& neiCc,
617 const labelList& testFaces,
618
619 labelList& globalRegion1,
620 labelList& globalRegion2
621 ) const;
622
623 //- Calculate intersections on zoned faces. Return per face -1
624 // or the global region of the surface and the orientation
625 // w.r.t. surface
626 void getIntersections
627 (
628 const labelList& surfacesToTest,
629 const pointField& neiCc,
630 const labelList& testFaces,
631
632 labelList& namedSurfaceRegion,
633 bitSet& posOrientation
634 ) const;
635
636 //- Determine patches for baffles
637 void getBafflePatches
638 (
639 const label nErodeCellZones,
640 const labelList& globalToMasterPatch,
641 const pointField& locationsInMesh,
642 const wordList& regionsInMesh,
643 const pointField& locationsOutsideMesh,
644 const bool exitIfLeakPath,
645 const refPtr<coordSetWriter>& leakPathFormatter,
646 refPtr<surfaceWriter>& surfFormatter,
647
648 const labelList& neiLevel,
649 const pointField& neiCc,
650 labelList& ownPatch,
651 labelList& neiPatch
652 ) const;
653
654 autoPtr<mapPolyMesh> splitMesh
655 (
656 const label nBufferLayers,
657 const labelList& globalToMasterPatch,
658 const labelList& globalToSlavePatch,
659 labelList& cellRegion,
660 labelList& ownPatch,
661 labelList& neiPatch
662 );
663
664 //- Repatches external face or creates baffle for internal face
665 // with user specified patches (might be different for both sides).
666 // Returns label of added face.
667 label createBaffle
668 (
669 const label facei,
670 const label ownPatch,
671 const label neiPatch,
672 polyTopoChange& meshMod
673 ) const;
674
675 //- Write leak path
676 static fileName writeLeakPath
677 (
678 const polyMesh& mesh,
679 const pointField& locationsInMesh,
680 const pointField& locationsOutsideMesh,
681 const boolList& blockedFace,
682 coordSetWriter& leakPathWriter
683 );
684
685
686 // Problem cell handling
687
688 //- Helper function to mark face as being on 'boundary'. Used by
689 // markFacesOnProblemCells
690 void markBoundaryFace
691 (
692 const label facei,
693 boolList& isBoundaryFace,
694 boolList& isBoundaryEdge,
695 boolList& isBoundaryPoint
696 ) const;
697
698 void findNearest
699 (
700 const labelList& meshFaces,
701 List<pointIndexHit>& nearestInfo,
702 labelList& nearestSurface,
703 labelList& nearestRegion,
704 vectorField& nearestNormal
705 ) const;
706
707 Map<label> findEdgeConnectedProblemCells
708 (
709 const scalarField& perpendicularAngle,
710 const labelList&
711 ) const;
712
713 bool isCollapsedFace
714 (
715 const pointField&,
716 const pointField& neiCc,
717 const scalar minFaceArea,
718 const scalar maxNonOrtho,
719 const label facei
720 ) const;
721
722 bool isCollapsedCell
723 (
724 const pointField&,
725 const scalar volFraction,
726 const label celli
727 ) const;
728
729 //- Returns list with for every internal face -1 or the patch
730 // they should be baffled into. If removeEdgeConnectedCells is set
731 // removes cells based on perpendicularAngle.
732 void markFacesOnProblemCells
733 (
734 const dictionary& motionDict,
735 const bool removeEdgeConnectedCells,
736 const scalarField& perpendicularAngle,
737 const labelList& globalToMasterPatch,
738
739 labelList& facePatch,
741 ) const;
742
743 //- Calculates for every face the nearest 'start' face. Any
744 // unreached face does not get set (faceToStart[facei] = -1)
745 void nearestFace
746 (
747 const labelUList& startFaces,
748 const bitSet& isBlockedFace,
749
751 labelList& faceToStart,
752 const label nIter = labelMax
753 ) const;
754
755 //- Calculates for every face the label of the nearest
756 // patch its zone. Any unreached face (disconnected mesh?) becomes
757 // adaptPatchIDs[0]
758 void nearestPatch
759 (
760 const labelList& adaptPatchIDs,
761 labelList& nearestPatch,
762 labelList& nearestZone
763 ) const;
764
765 //- Returns list with for every face the label of the nearest
766 // patch. Any unreached face (disconnected mesh?) becomes
767 // adaptPatchIDs[0]
768 labelList nearestPatch(const labelList& adaptPatchIDs) const;
769
770 //- Returns list with for every face the label of the nearest
771 // (global) region. Any unreached face (disconnected mesh?) becomes
772 // defaultRegion
773 labelList nearestIntersection
774 (
775 const labelList& surfacesToTest,
776 const label defaultRegion
777 ) const;
778
779 //- Returns list with for every internal face -1 or the patch
780 // they should be baffled into.
781 void markFacesOnProblemCellsGeometric
782 (
783 const snapParameters& snapParams,
784 const dictionary& motionDict,
785 const labelList& globalToMasterPatch,
786 const labelList& globalToSlavePatch,
787
788 labelList& facePatch,
790 ) const;
791
792
793 // Baffle merging
794
795 //- Extract those baffles (duplicate) faces that are on the edge
796 // of a baffle region. These are candidates for merging.
797 // samePatch : include only if both sides on the same patch.
798 List<labelPair> freeStandingBaffles
799 (
800 const List<labelPair>&,
801 const scalar freeStandingAngle,
802 const bool samePatch
803 ) const;
804
805
806 // Zone handling
807
808 //- Finds zone per cell for cells inside closed named surfaces.
809 // (uses geometric test for insideness)
810 // Adapts namedSurfaceRegion so all faces on boundary of cellZone
811 // have corresponding faceZone.
812 void findCellZoneGeometric
813 (
814 const pointField& neiCc,
815 const labelList& closedNamedSurfaces,
816 labelList& namedSurfaceRegion,
817 const labelList& surfaceToCellZone,
818 labelList& cellToZone
819 ) const;
820
821 //- Finds zone per cell for cells inside region for which name
822 // is specified.
823 void findCellZoneInsideWalk
824 (
825 const pointField& locationsInMesh,
826 const labelList& zonesInMesh,
827 const labelList& blockedFace, // per face -1 or some index >= 0
828 labelList& cellToZone
829 ) const;
830
831 //- Finds zone per cell for cells inside region for which name
832 // is specified.
833 void findCellZoneInsideWalk
834 (
835 const pointField& locationsInMesh,
836 const wordList& zoneNamesInMesh,
837 const labelList& faceToZone, // per face -1 or some index >= 0
838 labelList& cellToZone
839 ) const;
840
841 //- Determines cell zone from cell region information.
842 bool calcRegionToZone
843 (
844 const label backgroundZoneID,
845 const label surfZoneI,
846 const label ownRegion,
847 const label neiRegion,
848
849 labelList& regionToCellZone
850 ) const;
851
852 //- Finds zone per cell. Uses topological walk with all faces
853 // marked in unnamedSurfaceRegion (intersections with unnamed
854 // surfaces) and namedSurfaceRegion (intersections with named
855 // surfaces) regarded as blocked.
856 void findCellZoneTopo
857 (
858 const label backgroundZoneID,
859 const pointField& locationsInMesh,
860 const labelList& unnamedSurfaceRegion,
861 const labelList& namedSurfaceRegion,
862 const labelList& surfaceToCellZone,
863 labelList& cellToZone
864 ) const;
865
866 //- Opposite of findCellTopo: finds assigned cell connected to
867 // an unassigned one and puts it in the background zone.
868 void erodeCellZone
869 (
870 const label nErodeCellZones,
871 const label backgroundZoneID,
872 const labelList& unnamedSurfaceRegion,
873 const labelList& namedSurfaceRegion,
874 labelList& cellToZone
875 ) const;
876
877 //- Variation of findCellZoneTopo: walks from cellZones but ignores
878 // face intersections (unnamedSurfaceRegion). WIP
879 void growCellZone
880 (
881 const label nGrowCellZones,
882 const label backgroundZoneID,
883 labelList& unnamedSurfaceRegion1,
884 labelList& unnamedSurfaceRegion2,
885 labelList& namedSurfaceRegion,
886 labelList& cellToZone
887 ) const;
888
889 //- Make namedSurfaceRegion consistent with cellToZone
890 // - clear out any blocked faces inbetween same cell zone.
891 void makeConsistentFaceIndex
892 (
893 const labelList& zoneToNamedSurface,
894 const labelList& cellToZone,
895 labelList& namedSurfaceRegion
896 ) const;
897
898 //- Calculate cellZone allocation
899 void zonify
900 (
901 const bool allowFreeStandingZoneFaces,
902 const label nErodeCellZones,
903 const label backgroundZoneID,
904 const pointField& locationsInMesh,
905 const wordList& zonesInMesh,
906 const pointField& locationsOutsideMesh,
907 const bool exitIfLeakPath,
908 const refPtr<coordSetWriter>& leakPathFormatter,
909 refPtr<surfaceWriter>& surfFormatter,
910
911 labelList& cellToZone,
912 labelList& unnamedRegion1,
913 labelList& unnamedRegion2,
914 labelList& namedSurfaceRegion,
915 bitSet& posOrientation
916 ) const;
917
918 //- Put cells into cellZone, faces into faceZone
919 void zonify
920 (
921 const bitSet& isMasterFace,
922 const labelList& cellToZone,
923 const labelList& neiCellZone,
924 const labelList& faceToZone,
925 const bitSet& meshFlipMap,
926 polyTopoChange& meshMod
927 ) const;
928
929 //- Allocate faceZoneName
930 void allocateInterRegionFaceZone
931 (
932 const label ownZone,
933 const label neiZone,
934 wordPairHashTable& zonesToFaceZone,
935 LabelPairMap<word>& zoneIDsToFaceZone
936 ) const;
937
938 //- Remove any loose standing cells
939 void handleSnapProblems
940 (
941 const snapParameters& snapParams,
942 const bool useTopologicalSnapDetection,
943 const bool removeEdgeConnectedCells,
944 const scalarField& perpendicularAngle,
945 const dictionary& motionDict,
946 Time& runTime,
947 const labelList& globalToMasterPatch,
948 const labelList& globalToSlavePatch
949 );
950
951
952 // Some patch utilities
953
954 //- Get all faces in faceToZone that have no cellZone on
955 // either side.
956 labelList freeStandingBaffleFaces
957 (
958 const labelList& faceToZone,
959 const labelList& cellToZone,
960 const labelList& neiCellZone
961 ) const;
962
963 //- Determine per patch edge the number of master faces. Used
964 // to detect non-manifold situations.
965 void calcPatchNumMasterFaces
966 (
967 const bitSet& isMasterFace,
968 const indirectPrimitivePatch& patch,
969 labelList& nMasterFaces
970 ) const;
971
972 //- Determine per patch face the (singly-) connected zone it
973 // is in. Return overall number of zones.
974 label markPatchZones
975 (
976 const indirectPrimitivePatch& patch,
977 const labelList& nMasterFaces,
978 labelList& faceToZone
979 ) const;
980
981 //- Make faces consistent.
982 void consistentOrientation
983 (
984 const bitSet& isMasterFace,
985 const indirectPrimitivePatch& patch,
986 const labelList& nMasterFaces,
987 const labelList& faceToZone,
988 const Map<label>& zoneToOrientation,
989 bitSet& meshFlipMap
990 ) const;
991
992
993 //- No copy construct
994 meshRefinement(const meshRefinement&) = delete;
995
996 //- No copy assignment
997 void operator=(const meshRefinement&) = delete;
998
999public:
1000
1001 //- Runtime type information
1002 ClassName("meshRefinement");
1003
1004
1005 // Constructors
1006
1007 //- Construct from components
1008 meshRefinement
1009 (
1010 fvMesh& mesh,
1011 const scalar mergeDistance,
1012 const bool overwrite,
1013 const refinementSurfaces&,
1014 const refinementFeatures&,
1015 const shellSurfaces&, // omnidirectional refinement
1016 const shellSurfaces&, // limit refinement
1017 const labelUList& checkFaces, // initial faces to check
1018 const MeshType meshType,
1019 const bool dryRun
1020 );
1021
1022
1023 // Member Functions
1024
1025 // Access
1026
1027 //- Reference to mesh
1028 const fvMesh& mesh() const
1029 {
1030 return mesh_;
1031 }
1032 fvMesh& mesh()
1033 {
1034 return mesh_;
1035 }
1036
1037 scalar mergeDistance() const
1038 {
1039 return mergeDistance_;
1040 }
1041
1042 //- Overwrite the mesh?
1043 bool overwrite() const
1044 {
1045 return overwrite_;
1046 }
1047
1048 //- (points)instance of mesh upon construction
1049 const word& oldInstance() const
1050 {
1051 return oldInstance_;
1052 }
1053
1054 //- Reference to surface search engines
1055 const refinementSurfaces& surfaces() const
1056 {
1057 return surfaces_;
1058 }
1059
1060 //- Reference to feature edge mesh
1061 const refinementFeatures& features() const
1062 {
1063 return features_;
1064 }
1065
1066 //- Reference to refinement shells (regions)
1067 const shellSurfaces& shells() const
1068 {
1069 return shells_;
1070 }
1071
1072 //- Reference to limit shells (regions)
1073 const shellSurfaces& limitShells() const
1074 {
1075 return limitShells_;
1076 }
1077
1078 //- Reference to meshcutting engine
1079 const hexRef8& meshCutter() const
1080 {
1081 return meshCutter_;
1082 }
1083
1084 //- Mode of meshing
1085 MeshType meshType() const
1086 {
1087 return meshType_;
1088 }
1089
1090 //- Per start-end edge the index of the surface hit
1091 const labelList& surfaceIndex() const;
1092
1094
1095 //- For faces originating from processor faces store the original
1096 // patch
1097 const Map<label>& faceToCoupledPatch() const
1098 {
1099 return faceToCoupledPatch_;
1100 }
1101
1102 //- Additional face data that is maintained across
1103 // topo changes. Every entry is a list over all faces.
1104 // Bit of a hack. Additional flag to say whether to maintain master
1105 // only (false) or increase set to account for face-from-face.
1106 const List<Tuple2<mapType, labelList>>& userFaceData() const
1107 {
1108 return userFaceData_;
1109 }
1110
1111 List<Tuple2<mapType, labelList>>& userFaceData()
1112 {
1113 return userFaceData_;
1114 }
1115
1116
1117 // Other
1118
1119 //- Count number of intersections (local)
1120 label countHits() const;
1121
1122 //- Redecompose according to cell count
1123 // keepZoneFaces : find all faceZones from zoned surfaces and keep
1124 // owner and neighbour together
1125 // keepBaffles : find all baffles and keep them together
1126 // singleProcPoints : all faces using point are kept together
1127 autoPtr<mapDistributePolyMesh> balance
1128 (
1129 const bool keepZoneFaces,
1130 const bool keepBaffles,
1131 const labelList& singleProcPoints,
1132 const scalarField& cellWeights,
1133 decompositionMethod& decomposer,
1134 fvMeshDistribute& distributor
1135 );
1136
1137 //- Get faces with intersection.
1139
1140 //- Get points on surfaces with intersection and boundary faces.
1142
1143 //- Create patch from set of patches
1144 static autoPtr<indirectPrimitivePatch> makePatch
1145 (
1146 const polyMesh&,
1147 const labelList&
1148 );
1149
1150 //- Helper function to make a pointVectorField with correct
1151 // bcs for mesh movement:
1152 // - adaptPatchIDs : fixedValue
1153 // - processor : calculated (so free to move)
1154 // - cyclic/wedge/symmetry : slip
1155 // - other : slip
1156 static tmp<pointVectorField> makeDisplacementField
1157 (
1158 const pointMesh& pMesh,
1159 const labelList& adaptPatchIDs
1160 );
1161
1162 //- Helper function: check that face zones are synced
1163 static void checkCoupledFaceZones(const polyMesh&);
1164
1165 //- Helper: calculate edge weights (1/length)
1166 static void calculateEdgeWeights
1167 (
1168 const polyMesh& mesh,
1169 const bitSet& isMasterEdge,
1170 const labelList& meshPoints,
1171 const edgeList& edges,
1172 scalarField& edgeWeights,
1173 scalarField& invSumWeight
1174 );
1175
1176 //- Helper: weighted sum (over all subset of mesh points) by
1177 // summing contribution from (master) edges
1178 template<class Type>
1179 static void weightedSum
1180 (
1181 const polyMesh& mesh,
1182 const bitSet& isMasterEdge,
1183 const labelList& meshPoints,
1184 const edgeList& edges,
1185 const scalarField& edgeWeights,
1186 const Field<Type>& data,
1187 Field<Type>& sum
1188 );
1189
1190
1191 // Refinement
1192
1193 //- Is local topology a small gap?
1194 bool isGap
1195 (
1196 const scalar,
1197 const vector&,
1198 const vector&,
1199 const vector&,
1200 const vector&
1201 ) const;
1202
1203 //- Is local topology a small gap normal to the test vector
1204 bool isNormalGap
1205 (
1206 const scalar planarCos,
1207 const label level0,
1208 const vector& point0,
1209 const vector& normal0,
1210 const label level1,
1211 const vector& point1,
1212 const vector& normal1
1213 ) const;
1214
1215 //- Calculate list of cells to refine.
1217 (
1218 const pointField& keepPoints,
1219 const scalar curvature,
1220 const scalar planarAngle,
1221
1222 const bool featureRefinement,
1223 const bool featureDistanceRefinement,
1224 const bool internalRefinement,
1225 const bool surfaceRefinement,
1226 const bool curvatureRefinement,
1227 const bool smallFeatureRefinement,
1228 const bool gapRefinement,
1229 const bool bigGapRefinement,
1230 const bool spreadGapSize,
1231 const label maxGlobalCells,
1232 const label maxLocalCells
1233 ) const;
1234
1235
1236 // Blocking cells
1237
1238 //- Mark faces on interface between set and rest
1239 // (and same cell level)
1240 void markOutsideFaces
1241 (
1242 const labelList& cellLevel,
1243 const labelList& neiLevel,
1244 const labelList& refineCell,
1245 bitSet& isOutsideFace
1246 ) const;
1247
1248 //- Count number of faces on cell that are in set
1250 (
1251 const bitSet& isOutsideFace,
1252 const label celli
1253 ) const;
1255 //- Add one layer of cells to set
1256 void growSet
1257 (
1258 const labelList& neiLevel,
1259 const bitSet& isOutsideFace,
1261 label& nRefine
1262 ) const;
1263
1264 //- Detect gapRefinement cells and remove them
1266 (
1267 const scalar planarAngle,
1268 const labelList& minSurfaceLevel,
1269 const labelList& globalToMasterPatch,
1270 const label growIter
1271 );
1272
1273 // Blocking leaks (by blocking cells)
1274
1275 //- Faces currently on boundary or intersected by surface
1277 (
1280 ) const;
1281
1282 //- Return list of cells to block by walking from the seedCells
1283 // until reaching a leak face
1285 (
1287 const labelList& leakFaces,
1288 const labelList& seedCells
1289 ) const;
1290
1291 //- Remove minimum amount of cells to break any leak from
1292 // inside to outside
1295 const labelList& globalToMasterPatch,
1296 const labelList& globalToSlavePatch,
1297 const pointField& locationsInMesh,
1298 const pointField& locationsOutsideMesh,
1299 const labelList& selectedSurfaces
1300 );
1301
1302 //- Baffle faces to break any leak from inside to outside
1304 (
1305 const labelList& globalToMasterPatch,
1306 const labelList& globalToSlavePatch,
1307 const pointField& locationsInMesh,
1308 const wordList& zonesInMesh,
1309 const pointField& locationsOutsideMesh,
1310 const labelList& selectedSurfaces
1311 );
1312
1313
1314 //- Refine some cells
1315 autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
1316
1317 //- Balance the mesh
1319 (
1320 const string& msg,
1321 decompositionMethod& decomposer,
1322 fvMeshDistribute& distributor,
1323 labelList& cellsToRefine,
1324 const scalar maxLoadUnbalance,
1325 const label maxCellUnbalance
1326 );
1327
1328 //- Refine some cells and rebalance
1330 (
1331 const string& msg,
1332 decompositionMethod& decomposer,
1333 fvMeshDistribute& distributor,
1334 const labelList& cellsToRefine,
1335 const scalar maxLoadUnbalance,
1336 const label maxCellUnbalance
1337 );
1338
1339 //- Balance before refining some cells
1341 (
1342 const string& msg,
1343 decompositionMethod& decomposer,
1344 fvMeshDistribute& distributor,
1345 const labelList& cellsToRefine,
1346 const scalar maxLoadUnbalance,
1347 const label maxCellUnbalance
1348 );
1349
1350 //- Calculate list of cells to directionally refine
1353 const label maxGlobalCells,
1354 const label maxLocalCells,
1355 const labelList& currentLevel,
1356 const direction dir
1357 ) const;
1358
1359 //- Directionally refine in direction cmpt
1361 (
1362 const string& msg,
1363 const direction cmpt,
1364 const labelList& cellsToRefine
1365 );
1366
1367
1368 // Baffle handling
1369
1370 //- Split off unreachable areas of mesh.
1372 (
1373 const bool handleSnapProblems,
1374
1375 // How to remove problem snaps
1376 const snapParameters& snapParams,
1377 const bool useTopologicalSnapDetection,
1378 const bool removeEdgeConnectedCells,
1379 const scalarField& perpendicularAngle,
1380 const label nErodeCellZones,
1381 const dictionary& motionDict,
1382 Time& runTime,
1383 const labelList& globalToMasterPatch,
1384 const labelList& globalToSlavePatch,
1385 const pointField& locationsInMesh,
1386 const wordList& regionsInMesh,
1387 const pointField& locationsOutsideMesh,
1388 const bool exitIfLeakPath,
1389 const refPtr<coordSetWriter>& leakPathFormatter,
1390 refPtr<surfaceWriter>& surfFormatter
1391 );
1392
1393 //- Merge free-standing baffles
1395 (
1396 const bool samePatch,
1397 const snapParameters& snapParams,
1398 const bool useTopologicalSnapDetection,
1399 const bool removeEdgeConnectedCells,
1400 const scalarField& perpendicularAngle,
1401 const scalar planarAngle,
1402 const dictionary& motionDict,
1403 Time& runTime,
1404 const labelList& globalToMasterPatch,
1405 const labelList& globalToSlavePatch,
1406 const pointField& locationsInMesh,
1407 const pointField& locationsOutsideMesh
1408 );
1409
1410 //- Split off (with optional buffer layers) unreachable areas
1411 // of mesh. Does not introduce baffles.
1412 autoPtr<mapPolyMesh> splitMesh
1413 (
1414 const label nBufferLayers,
1415 const label nErodeCellZones,
1416 const labelList& globalToMasterPatch,
1417 const labelList& globalToSlavePatch,
1418
1419 const pointField& locationsInMesh,
1420 const wordList& regionsInMesh,
1421 const pointField& locationsOutsideMesh,
1422 const bool exitIfLeakPath,
1423 const refPtr<coordSetWriter>& leakPathFormatter,
1424 refPtr<surfaceWriter>& surfFormatter
1425 );
1426
1427 //- Remove cells from limitRegions if level -1
1429 (
1430 const label nBufferLayers,
1431 const label nErodeCellZones,
1432 const labelList& globalToMasterPatch,
1433 const labelList& globalToSlavePatch,
1434 const pointField& locationsInMesh,
1435 const wordList& regionsInMesh,
1436 const pointField& locationsOutsideMesh
1437 );
1438
1439 //- Find boundary points that connect to more than one cell
1440 // region and split them.
1442
1443 //- Find boundary points that connect to more than one cell
1444 // region and split them.
1446
1447 //- Find boundary points that are on faceZones of type boundary
1448 // and duplicate them
1450
1451 //- Merge duplicate points
1453 (
1454 const labelList& pointToDuplicate
1455 );
1456
1457 //- Create baffle for every internal face where ownPatch != -1.
1458 // External faces get repatched according to ownPatch (neiPatch
1459 // should be -1 for these)
1461 (
1462 const labelList& ownPatch,
1463 const labelList& neiPatch
1464 );
1465
1466 //- Get zones of given type
1468 (
1470 ) const;
1471
1472 //- Subset baffles according to zones
1474 (
1475 const polyMesh& mesh,
1476 const labelList& zoneIDs,
1477 const List<labelPair>& baffles
1478 );
1479
1480 //- Map baffles after layer addition. Gets new-to-old face map.
1481 static void mapBaffles
1482 (
1483 const polyMesh& mesh,
1484 const labelList& faceMap,
1485 List<labelPair>& baffles
1486 );
1487
1488 //- Get per-face information (faceZone, master/slave patch)
1489 void getZoneFaces
1490 (
1491 const labelList& zoneIDs,
1493 labelList& ownPatch,
1494 labelList& neiPatch,
1495 labelList& nBaffles
1496 ) const;
1497
1498 //- Create baffles for faces on faceZones. Return created baffles
1499 // (= pairs of faces) and corresponding faceZone
1501 (
1502 const labelList& zoneIDs,
1503 List<labelPair>& baffles,
1504 labelList& originatingFaceZone
1505 );
1506
1507 //- Merge baffles. Gets pairs of faces and boundary faces to move
1508 // onto (coupled) patches
1510 (
1511 const List<labelPair>&,
1512 const Map<label>& faceToPatch
1513 );
1514
1515 //- Merge all baffles on faceZones
1517 (
1518 const bool doInternalZones,
1519 const bool doBaffleZones
1520 );
1521
1522 //- Put faces/cells into zones according to surface specification.
1523 // Returns null if no zone surfaces present. Regions containing
1524 // locationsInMesh/regionsInMesh will be put in corresponding
1525 // cellZone. keepPoints is for backwards compatibility and sets
1526 // all yet unassigned cells to be non-zoned (zone = -1)
1528 (
1529 const bool allowFreeStandingZoneFaces,
1530 const label nErodeCellZones,
1531 const pointField& locationsInMesh,
1532 const wordList& regionsInMesh,
1533 const pointField& locationsOutsideMesh,
1534 const bool exitIfLeakPath,
1535 const refPtr<coordSetWriter>& leakPathFormatter,
1536 refPtr<surfaceWriter>& surfFormatter,
1537 wordPairHashTable& zonesToFaceZone
1538 );
1539
1540
1541 // Other topo changes
1542
1543 //- Helper:append patch to end of mesh.
1544 static label appendPatch
1545 (
1546 fvMesh&,
1547 const label insertPatchi,
1548 const word&,
1549 const dictionary&
1550 );
1551
1552 //- Helper:add patch to mesh. Update all registered fields.
1553 // Used by addMeshedPatch to add patches originating from surfaces.
1554 static label addPatch(fvMesh&, const word& name, const dictionary&);
1555
1556 //- Add patch originating from meshing. Update meshedPatches_.
1557 label addMeshedPatch(const word& name, const dictionary&);
1558
1559 //- Get patchIDs for patches added in addMeshedPatch.
1560 labelList meshedPatches() const;
1561
1562 //- Add/lookup faceZone and update information. Return index of
1563 // faceZone
1564 label addFaceZone
1565 (
1566 const word& fzName,
1567 const word& masterPatch,
1568 const word& slavePatch,
1569 const surfaceZonesInfo::faceZoneType& fzType
1570 );
1571
1572 //- Lookup faceZone information. Return false if no information
1573 // for faceZone
1574 bool getFaceZoneInfo
1575 (
1576 const word& fzName,
1577 label& masterPatchID,
1578 label& slavePatchID,
1580 ) const;
1581
1582 //- Add pointZone if does not exist. Return index of zone
1583 label addPointZone(const word& name);
1584
1585 //- Count number of faces per patch edge. Parallel consistent.
1587
1588 //- Select coupled faces that are not collocated
1590
1591 //- Find any intersection of surface. Store in surfaceIndex_.
1592 void updateIntersections(const labelUList& changedFaces);
1593
1594 //- Calculate nearest intersection for selected mesh faces
1595 void nearestIntersection
1596 (
1597 const labelList& surfacesToTest,
1598 const labelList& testFaces,
1599
1600 labelList& surface1,
1601 List<pointIndexHit>& hit1,
1602 labelList& region1,
1603 labelList& surface2,
1604 List<pointIndexHit>& hit2,
1605 labelList& region2
1606 ) const;
1607
1608 //- Remove cells. Put exposedFaces into exposedPatchIDs.
1610 (
1611 const labelList& cellsToRemove,
1612 const labelList& exposedFaces,
1613 const labelList& exposedPatchIDs,
1614 removeCells& cellRemover
1615 );
1616
1617 //- Find cell point is in. Uses optional perturbation to re-test.
1618 // Returns -1 on processors that do not have the cell.
1619 static label findCell
1620 (
1621 const polyMesh&,
1622 const vector& perturbVec,
1623 const point& p
1624 );
1625
1626 //- Find region point is in. Uses optional perturbation to re-test.
1627 static label findRegion
1628 (
1629 const polyMesh&,
1630 const labelList& cellRegion,
1631 const vector& perturbVec,
1632 const point& p
1633 );
1634
1635 //- Find regions points are in.
1636 // \return number of cells to be removed
1637 static label findRegions
1638 (
1639 const polyMesh&,
1640 const vector& perturbVec,
1641 const pointField& locationsInMesh,
1642 const pointField& locationsOutsideMesh,
1643 const label nRegions,
1644 labelList& cellRegion,
1645 const boolList& blockedFace,
1646 // Leak-path
1647 const bool exitIfLeakPath,
1648 const refPtr<coordSetWriter>& leakPathFormatter
1649 );
1650
1651 //- Split mesh. Keep part containing point. Return empty map if
1652 // no cells removed.
1654 (
1655 const labelList& globalToMasterPatch,
1656 const labelList& globalToSlavePatch,
1657 const pointField& locationsInMesh,
1658 const pointField& locationsOutsideMesh,
1659 // Leak-path
1660 const bool exitIfLeakPath,
1661 const refPtr<coordSetWriter>& leakPathFormatter
1662 );
1663
1664 //- Helper: split face into:
1665 // - f0 : split[0] to split[1]
1666 // - f1 : split[1] to split[0]
1667 static void splitFace
1668 (
1669 const face& f,
1670 const labelPair& split,
1671
1672 face& f0,
1673 face& f1
1674 );
1675
1676 //- Split faces into two
1677 void doSplitFaces
1678 (
1679 const labelList& splitFaces,
1680 const labelPairList& splits,
1681 const labelPairList& splitPatches,
1682 polyTopoChange& meshMod
1683 ) const;
1684
1685 //- Split faces along diagonal. Maintain mesh quality. Return
1686 // total number of faces split.
1687 label splitFacesUndo
1688 (
1689 const labelList& splitFaces,
1690 const labelPairList& splits,
1691 const labelPairList& splitPatches,
1692 const dictionary& motionDict,
1693
1694 labelList& duplicateFace,
1695 List<labelPair>& baffles
1696 );
1697
1698 //- Update local numbering for mesh redistribution
1699 void distribute(const mapDistributePolyMesh&);
1700
1701 //- Update for external change to mesh. changedFaces are in new mesh
1702 // face labels.
1703 void updateMesh
1704 (
1705 const mapPolyMesh&,
1706 const labelList& changedFaces
1707 );
1708
1709 //- Helper: reorder list according to map.
1710 template<class T>
1711 static void updateList
1712 (
1713 const labelList& newToOld,
1714 const T& nullValue,
1715 List<T>& elems
1716 );
1717
1718
1719 // Restoring : is where other processes delete and reinsert data.
1720
1721 //- Signal points/face/cells for which to store data
1722 void storeData
1723 (
1724 const labelList& pointsToStore,
1725 const labelList& facesToStore,
1726 const labelList& cellsToStore
1727 );
1728
1729 //- Update local numbering + undo
1730 // Data to restore given as new pointlabel + stored pointlabel
1731 // (i.e. what was in pointsToStore)
1732 void updateMesh
1733 (
1734 const mapPolyMesh&,
1735 const labelList& changedFaces,
1736 const Map<label>& pointsToRestore,
1737 const Map<label>& facesToRestore,
1738 const Map<label>& cellsToRestore
1739 );
1740
1741 // Merging coplanar faces and edges
1742
1743 //- Merge coplanar faces if sets are of size mergeSize
1744 // (usually 4)
1745 label mergePatchFaces
1746 (
1747 const scalar minCos,
1748 const scalar concaveCos,
1749 const label mergeSize,
1750 const labelList& patchIDs,
1751 const meshRefinement::FaceMergeType mergeType
1752 );
1753
1754 //- Merge coplanar faces. preserveFaces is != -1 for faces
1755 // to be preserved
1757 (
1758 const scalar minCos,
1759 const scalar concaveCos,
1760 const labelList& patchIDs,
1761 const dictionary& motionDict,
1762 const labelList& preserveFaces,
1763 const meshRefinement::FaceMergeType mergeType
1764 );
1765
1767 (
1768 removePoints& pointRemover,
1769 const boolList& pointCanBeDeleted
1770 );
1771
1773 (
1774 removePoints& pointRemover,
1775 const labelList& facesToRestore
1776 );
1777
1779 (
1780 const labelList& candidateFaces,
1781 const labelHashSet& set
1782 ) const;
1783
1784 // Pick up faces of cells of faces in set.
1786 (
1787 const labelUList& set
1788 ) const;
1789
1790 // Pick up faces of cells of faces in set.
1792 (
1793 const labelHashSet& set
1794 ) const;
1795
1796 //- Merge edges, maintain mesh quality. Return global number
1797 // of edges merged
1798 label mergeEdgesUndo
1799 (
1800 const scalar minCos,
1801 const dictionary& motionDict
1802 );
1803
1804
1805 // Debug/IO
1806
1807 //- Debugging: check that all faces still obey start()>end()
1808 void checkData();
1809
1810 static void testSyncPointList
1811 (
1812 const string& msg,
1813 const polyMesh& mesh,
1814 const List<scalar>& fld
1815 );
1816
1817 static void testSyncPointList
1818 (
1819 const string& msg,
1820 const polyMesh& mesh,
1821 const List<point>& fld
1822 );
1823
1824 //- Compare two lists over all boundary faces
1825 template<class T>
1827 (
1828 const scalar mergeDistance,
1829 const string&,
1830 const UList<T>&,
1831 const UList<T>&
1832 ) const;
1833
1834 //- Print list according to (collected and) sorted coordinate
1835 template<class T>
1836 static void collectAndPrint
1837 (
1838 const UList<point>& points,
1839 const UList<T>& data
1840 );
1841
1842 //- Determine master point for subset of points. If coupled
1843 // chooses only one
1844 static bitSet getMasterPoints
1845 (
1846 const polyMesh& mesh,
1847 const labelList& meshPoints
1848 );
1849
1850 //- Determine master edge for subset of edges. If coupled
1851 // chooses only one
1852 static bitSet getMasterEdges
1853 (
1854 const polyMesh& mesh,
1855 const labelList& meshEdges
1856 );
1857
1858 //- Print some mesh stats.
1859 void printMeshInfo
1860 (
1861 const bool debug,
1862 const string& msg,
1863 const bool printCellLevel
1864 ) const;
1865
1866 //- Replacement for Time::timeName() that returns oldInstance
1867 //- (if overwrite_)
1868 word timeName() const;
1869
1870 //- Set instance of all local IOobjects
1871 void setInstance(const fileName&);
1872
1873 //- Write mesh and all data
1874 bool write() const;
1875
1876 //- Write refinement level as volScalarFields for postprocessing
1877 void dumpRefinementLevel() const;
1878
1879 //- Debug: Write intersection information to OBJ format
1880 void dumpIntersections(const fileName& prefix) const;
1881
1882 //- Do any one of above IO functions
1883 void write
1884 (
1885 const debugType debugFlags,
1886 const writeType writeFlags,
1887 const fileName&
1888 ) const;
1889
1890 //- Helper: remove all relevant files from mesh instance
1891 static void removeFiles(const polyMesh&);
1892
1893 //- Helper: calculate average
1894 template<class T>
1895 static T gAverage
1896 (
1897 const bitSet& isMasterElem,
1898 const UList<T>& values
1899 );
1900
1901 //- Get/set write level
1902 static writeType writeLevel();
1903 static void writeLevel(const writeType);
1904
1906 //static outputType outputLevel();
1907 //static void outputLevel(const outputType);
1908
1909
1910 //- Helper: convert wordList into bit pattern using provided Enum
1911 template<class EnumContainer>
1912 static int readFlags
1913 (
1914 const EnumContainer& namedEnum,
1915 const wordList& words
1916 );
1917
1918 //- Wrapper around dictionary::get which does not exit
1919 template<class Type>
1920 static Type get
1921 (
1922 const dictionary& dict,
1923 const word& keyword,
1924 const bool noExit,
1925 enum keyType::option matchOpt = keyType::REGEX,
1926 const Type& deflt = Zero
1927 );
1928
1929 //- Wrapper around dictionary::subDict which does not exit
1931 static const dictionary& subDict
1932 (
1933 const dictionary& dict,
1934 const word& keyword,
1935 const bool noExit,
1936 enum keyType::option matchOpt = keyType::REGEX
1937 );
1938
1939 //- Wrapper around dictionary::lookup which does not exit
1940 static ITstream& lookup
1941 (
1942 const dictionary& dict,
1943 const word& keyword,
1944 const bool noExit,
1945 enum keyType::option matchOpt = keyType::REGEX
1946 );
1947};
1948
1949
1950// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1951
1952} // End namespace Foam
1953
1954// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1955
1956#ifdef NoRepository
1957 #include "meshRefinementTemplates.C"
1958#endif
1959
1960// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1961
1962#endif
1963
1964// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
An input stream of tokens.
Definition ITstream.H:56
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 HashTable to objects of type <T> with a label key.
Definition Map.H:54
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
Base class for writing coordSet(s) and tracks with fields.
Abstract base class for domain decomposition.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
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
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Refinement of (split) hexes using polyTopoChange.
Definition hexRef8.H:63
option
Enumeration for the data type and search/match modes (bitmask).
Definition keyType.H:80
@ REGEX
Regular expression.
Definition keyType.H:83
Takes mesh with 'baffles' (= boundary faces sharing points). Determines for selected points on bounda...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing processor-to-processor mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
static void calculateEdgeWeights(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length).
label addFaceZone(const word &fzName, const word &masterPatch, const word &slavePatch, const surfaceZonesInfo::faceZoneType &fzType)
Add/lookup faceZone and update information. Return index of.
labelList getZones(const List< surfaceZonesInfo::faceZoneType > &fzTypes) const
Get zones of given type.
bool getFaceZoneInfo(const word &fzName, label &masterPatchID, label &slavePatchID, surfaceZonesInfo::faceZoneType &fzType) const
Lookup faceZone information. Return false if no information.
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
static void collectAndPrint(const UList< point > &points, const UList< T > &data)
Print list according to (collected and) sorted coordinate.
labelList intersectedFaces() const
Get faces with intersection.
labelList intersectedPoints() const
Get points on surfaces with intersection and boundary faces.
static ITstream & lookup(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Wrapper around dictionary::lookup which does not exit.
autoPtr< mapPolyMesh > removeLimitShells(const label nBufferLayers, const label nErodeCellZones, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList &regionsInMesh, const pointField &locationsOutsideMesh)
Remove cells from limitRegions if level -1.
void printMeshInfo(const bool debug, const string &msg, const bool printCellLevel) const
Print some mesh stats.
MeshType
Enumeration for how to operate.
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
label mergeEdgesUndo(const scalar minCos, const dictionary &motionDict)
Merge edges, maintain mesh quality. Return global number.
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
static void mapBaffles(const polyMesh &mesh, const labelList &faceMap, List< labelPair > &baffles)
Map baffles after layer addition. Gets new-to-old face map.
static const Enum< MeshType > MeshTypeNames
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
void checkData()
Debugging: check that all faces still obey start()>end().
autoPtr< mapPolyMesh > blockLeakFaces(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList &zonesInMesh, const pointField &locationsOutsideMesh, const labelList &selectedSurfaces)
Baffle faces to break any leak from inside to outside.
static List< labelPair > subsetBaffles(const polyMesh &mesh, const labelList &zoneIDs, const List< labelPair > &baffles)
Subset baffles according to zones.
const List< Tuple2< mapType, labelList > > & userFaceData() const
Additional face data that is maintained across.
void updateIntersections(const labelUList &changedFaces)
Find any intersection of surface. Store in surfaceIndex_.
static Type get(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX, const Type &deflt=Zero)
Wrapper around dictionary::get which does not exit.
void updateMesh(const mapPolyMesh &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
writeType
Enumeration for what to write. Used as a bit-pattern.
void doSplitFaces(const labelList &splitFaces, const labelPairList &splits, const labelPairList &splitPatches, polyTopoChange &meshMod) const
Split faces into two.
void mergeFreeStandingBaffles(const bool samePatch, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const scalar planarAngle, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh)
Merge free-standing baffles.
const hexRef8 & meshCutter() const
Reference to meshcutting engine.
autoPtr< mapPolyMesh > mergePoints(const labelList &pointToDuplicate)
Merge duplicate points.
const labelList & surfaceIndex() const
Per start-end edge the index of the surface hit.
autoPtr< mapPolyMesh > doRestorePoints(removePoints &pointRemover, const labelList &facesToRestore)
autoPtr< mapPolyMesh > dupNonManifoldPoints()
Find boundary points that connect to more than one cell.
labelList refineCandidates(const pointField &keepPoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool smallFeatureRefinement, const bool gapRefinement, const bool bigGapRefinement, const bool spreadGapSize, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
const shellSurfaces & limitShells() const
Reference to limit shells (regions).
label countFaceDirs(const bitSet &isOutsideFace, const label celli) const
Count number of faces on cell that are in set.
labelList growFaceCellFace(const labelUList &set) const
void dumpRefinementLevel() const
Write refinement level as volScalarFields for postprocessing.
autoPtr< mapPolyMesh > removeGapCells(const scalar planarAngle, const labelList &minSurfaceLevel, const labelList &globalToMasterPatch, const label growIter)
Detect gapRefinement cells and remove them.
autoPtr< mapPolyMesh > dupNonManifoldBoundaryPoints()
Find boundary points that are on faceZones of type boundary.
static const Enum< writeType > writeTypeNames
static int readFlags(const EnumContainer &namedEnum, const wordList &words)
Helper: convert wordList into bit pattern using provided Enum.
void distribute(const mapDistributePolyMesh &)
Update local numbering for mesh redistribution.
static void splitFace(const face &f, const labelPair &split, face &f0, face &f1)
Helper: split face into:
void getZoneFaces(const labelList &zoneIDs, labelList &faceZoneID, labelList &ownPatch, labelList &neiPatch, labelList &nBaffles) const
Get per-face information (faceZone, master/slave patch).
const Map< label > & faceToCoupledPatch() const
For faces originating from processor faces store the original.
autoPtr< mapPolyMesh > mergeBaffles(const List< labelPair > &, const Map< label > &faceToPatch)
Merge baffles. Gets pairs of faces and boundary faces to move.
static label addPatch(fvMesh &, const word &name, const dictionary &)
Helper:add patch to mesh. Update all registered fields.
label countHits() const
Count number of intersections (local).
void markOutsideFaces(const labelList &cellLevel, const labelList &neiLevel, const labelList &refineCell, bitSet &isOutsideFace) const
Mark faces on interface between set and rest.
ClassName("meshRefinement")
Runtime type information.
autoPtr< mapPolyMesh > doRemoveCells(const labelList &cellsToRemove, const labelList &exposedFaces, const labelList &exposedPatchIDs, removeCells &cellRemover)
Remove cells. Put exposedFaces into exposedPatchIDs.
mapType
Enumeration for how userdata is to be mapped upon refinement.
@ KEEPALL
have slaves (upon refinement) from master
@ REMOVE
set value to -1 any face that was refined
@ MASTERONLY
maintain master only
static FOAM_NO_DANGLING_REFERENCE const dictionary & subDict(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Wrapper around dictionary::subDict which does not exit.
static label findCell(const polyMesh &, const vector &perturbVec, const point &p)
Find cell point is in. Uses optional perturbation to re-test.
const refinementSurfaces & surfaces() const
Reference to surface search engines.
label mergePatchFaces(const scalar minCos, const scalar concaveCos, const label mergeSize, const labelList &patchIDs, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces if sets are of size mergeSize.
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_).
debugType
Enumeration for what to debug. Used as a bit-pattern.
static bitSet getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
const shellSurfaces & shells() const
Reference to refinement shells (regions).
labelList meshedPatches() const
Get patchIDs for patches added in addMeshedPatch.
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
labelList countEdgeFaces(const uindirectPrimitivePatch &pp) const
Count number of faces per patch edge. Parallel consistent.
bool isNormalGap(const scalar planarCos, const label level0, const vector &point0, const vector &normal0, const label level1, const vector &point1, const vector &normal1) const
Is local topology a small gap normal to the test vector.
labelList collectFaces(const labelList &candidateFaces, const labelHashSet &set) const
static void weightedSum(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
autoPtr< mapPolyMesh > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Split mesh. Keep part containing point. Return empty map if.
autoPtr< mapDistributePolyMesh > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Refine some cells and rebalance.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
List< Tuple2< mapType, labelList > > & userFaceData()
const fvMesh & mesh() const
Reference to mesh.
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
const refinementFeatures & features() const
Reference to feature edge mesh.
void testSyncBoundaryFaceList(const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
Compare two lists over all boundary faces.
void selectSeparatedCoupledFaces(boolList &) const
Select coupled faces that are not collocated.
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
const word & oldInstance() const
(points)instance of mesh upon construction
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &p)
Find region point is in. Uses optional perturbation to re-test.
label addPointZone(const word &name)
Add pointZone if does not exist. Return index of zone.
void selectIntersectedFaces(const labelList &surfaces, boolList &isBlockedFace) const
Faces currently on boundary or intersected by surface.
static label findRegions(const polyMesh &, const vector &perturbVec, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const label nRegions, labelList &cellRegion, const boolList &blockedFace, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Find regions points are in.
autoPtr< mapPolyMesh > doRemovePoints(removePoints &pointRemover, const boolList &pointCanBeDeleted)
scalar mergeDistance() const
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
autoPtr< mapPolyMesh > directionalRefine(const string &msg, const direction cmpt, const labelList &cellsToRefine)
Directionally refine in direction cmpt.
autoPtr< mapPolyMesh > refine(const labelList &cellsToRefine)
Refine some cells.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
void storeData(const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
Signal points/face/cells for which to store data.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
bool overwrite() const
Overwrite the mesh?
bool write() const
Write mesh and all data.
autoPtr< mapPolyMesh > mergeZoneBaffles(const bool doInternalZones, const bool doBaffleZones)
Merge all baffles on faceZones.
labelList directionalRefineCandidates(const label maxGlobalCells, const label maxLocalCells, const labelList &currentLevel, const direction dir) const
Calculate list of cells to directionally refine.
void growSet(const labelList &neiLevel, const bitSet &isOutsideFace, labelList &refineCell, label &nRefine) const
Add one layer of cells to set.
autoPtr< mapDistributePolyMesh > balance(const bool keepZoneFaces, const bool keepBaffles, const labelList &singleProcPoints, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
label mergePatchFacesUndo(const scalar minCos, const scalar concaveCos, const labelList &patchIDs, const dictionary &motionDict, const labelList &preserveFaces, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces. preserveFaces is != -1 for faces.
label splitFacesUndo(const labelList &splitFaces, const labelPairList &splits, const labelPairList &splitPatches, const dictionary &motionDict, labelList &duplicateFace, List< labelPair > &baffles)
Split faces along diagonal. Maintain mesh quality. Return.
autoPtr< mapPolyMesh > removeLeakCells(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const labelList &selectedSurfaces)
Remove minimum amount of cells to break any leak from.
void baffleAndSplitMesh(const bool handleSnapProblems, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const label nErodeCellZones, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList &regionsInMesh, const pointField &locationsOutsideMesh, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter, refPtr< surfaceWriter > &surfFormatter)
Split off unreachable areas of mesh.
void setInstance(const fileName &)
Set instance of all local IOobjects.
autoPtr< mapPolyMesh > createBaffles(const labelList &ownPatch, const labelList &neiPatch)
Create baffle for every internal face where ownPatch != -1.
static const Enum< debugType > debugTypeNames
autoPtr< mapDistributePolyMesh > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Balance before refining some cells.
static label appendPatch(fvMesh &, const label insertPatchi, const word &, const dictionary &)
Helper:append patch to end of mesh.
MeshType meshType() const
Mode of meshing.
static writeType writeLevel()
Get/set write level.
autoPtr< mapPolyMesh > createZoneBaffles(const labelList &zoneIDs, List< labelPair > &baffles, labelList &originatingFaceZone)
Create baffles for faces on faceZones. Return created baffles.
void dumpIntersections(const fileName &prefix) const
Debug: Write intersection information to OBJ format.
labelList detectLeakCells(const boolList &isBlockedFace, const labelList &leakFaces, const labelList &seedCells) const
Return list of cells to block by walking from the seedCells.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
Container with cells to refine. Refinement given as single direction.
Definition refineCell.H:53
Encapsulates queries for features.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
Given list of cells to remove, insert all the topology changes.
Definition removeCells.H:60
Removes selected points from mesh and updates faces using these points.
Encapsulates queries for volume refinement ('refine all cells within shell').
Simple container to keep together snap specific information.
Contains information about location on a triSurface.
faceZoneType
What to do with faceZone faces.
A class for managing temporary objects.
Definition tmp.H:75
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition volumeType.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
static bool split(const std::string &line, std::string &key, std::string &val)
Definition cpuInfo.C:32
dynamicFvMesh & mesh
engineTime & runTime
const pointField & points
const labelIOList & zoneIDs
Definition correctPhi.H:59
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
const bitSet isBlockedFace(intersectedFaces())
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
HashTable< word, wordPair, Foam::Hash< wordPair > > wordPairHashTable
HashTable of wordPair.
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
HashTable< T, labelPair, Foam::Hash< labelPair > > LabelPairMap
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
DynamicID< faceZoneMesh > faceZoneID
Foam::faceZoneID.
Definition ZoneIDs.H:43
Forwards and collection of common point field types.
labelList f(nPoints)
dictionary dict
volScalarField & e
#define FOAM_NO_DANGLING_REFERENCE
Definition stdFoam.H:80