Loading...
Searching...
No Matches
snappyRefineDriver.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-2015 OpenFOAM Foundation
9 Copyright (C) 2015-2025 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 "snappyRefineDriver.H"
30#include "meshRefinement.H"
31#include "fvMesh.H"
32#include "Time.H"
33#include "cellSet.H"
34#include "syncTools.H"
36#include "refinementSurfaces.H"
37#include "refinementFeatures.H"
38#include "shellSurfaces.H"
40#include "unitConversion.H"
41#include "snapParameters.H"
42#include "localPointRegion.H"
43#include "IOmanip.H"
44#include "labelVector.H"
45#include "profiling.H"
46#include "searchableSurfaces.H"
47#include "fvMeshSubset.H"
48#include "interpolationTable.H"
50#include "regionSplit.H"
51#include "removeCells.H"
52#include "addPatchCellLayer.H"
53
54// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
55
56namespace Foam
58 defineTypeNameAndDebug(snappyRefineDriver, 0);
59} // End namespace Foam
60
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
64Foam::snappyRefineDriver::snappyRefineDriver
65(
66 meshRefinement& meshRefiner,
67 decompositionMethod& decomposer,
68 fvMeshDistribute& distributor,
69 const labelUList& globalToMasterPatch,
70 const labelUList& globalToSlavePatch,
71 coordSetWriter& setFormatter,
72 refPtr<surfaceWriter>& surfFormatter,
73 const bool dryRun
74)
75:
76 meshRefiner_(meshRefiner),
77 decomposer_(decomposer),
78 distributor_(distributor),
79 globalToMasterPatch_(globalToMasterPatch),
80 globalToSlavePatch_(globalToSlavePatch),
81 setFormatter_(setFormatter),
82 surfFormatter_(surfFormatter),
83 dryRun_(dryRun)
84{}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
89Foam::label Foam::snappyRefineDriver::featureEdgeRefine
90(
91 const refinementParameters& refineParams,
92 const label maxIter,
93 const label minRefine
94)
95{
96 if (dryRun_)
97 {
98 return 0;
99 }
100
101 if (refineParams.minRefineCells() == -1)
102 {
103 // Special setting to be able to restart shm on meshes with inconsistent
104 // cellLevel/pointLevel
105 return 0;
106 }
107
108 addProfiling(edge, "snappyHexMesh::refine::edge");
109 const fvMesh& mesh = meshRefiner_.mesh();
110
111 label iter = 0;
112
113 if (meshRefiner_.features().size() && maxIter > 0)
114 {
115 for (; iter < maxIter; iter++)
116 {
117 Info<< nl
118 << "Feature refinement iteration " << iter << nl
119 << "------------------------------" << nl
120 << endl;
121
122 labelList candidateCells
123 (
124 meshRefiner_.refineCandidates
125 (
126 refineParams.locationsInMesh(),
127 refineParams.curvature(),
128 refineParams.planarAngle(),
129
130 true, // featureRefinement
131 false, // featureDistanceRefinement
132 false, // internalRefinement
133 false, // surfaceRefinement
134 false, // curvatureRefinement
135 false, // smallFeatureRefinement
136 false, // gapRefinement
137 false, // bigGapRefinement
138 false, // spreadGapSize
139 refineParams.maxGlobalCells(),
140 refineParams.maxLocalCells()
141 )
142 );
143 labelList cellsToRefine
144 (
145 meshRefiner_.meshCutter().consistentRefinement
146 (
147 candidateCells,
148 true
149 )
150 );
151 Info<< "Determined cells to refine in = "
152 << mesh.time().cpuTimeIncrement() << " s" << endl;
153
154
155
156 const label nCellsToRefine =
157 returnReduce(cellsToRefine.size(), sumOp<label>());
158
159 Info<< "Selected for feature refinement : " << nCellsToRefine
160 << " cells (out of " << mesh.globalData().nTotalCells()
161 << ')' << endl;
162
163 if (nCellsToRefine <= minRefine)
164 {
165 Info<< "Stopping refining since too few cells selected."
166 << nl << endl;
167 break;
168 }
169
170
171 if (debug > 0)
172 {
173 const_cast<Time&>(mesh.time())++;
174 }
175
176
177 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
178 {
179 meshRefiner_.balanceAndRefine
180 (
181 "feature refinement iteration " + name(iter),
182 decomposer_,
183 distributor_,
184 cellsToRefine,
185 refineParams.maxLoadUnbalance(),
186 refineParams.maxCellUnbalance()
187 );
188 }
189 else
190 {
191 meshRefiner_.refineAndBalance
192 (
193 "feature refinement iteration " + name(iter),
194 decomposer_,
195 distributor_,
196 cellsToRefine,
197 refineParams.maxLoadUnbalance(),
198 refineParams.maxCellUnbalance()
199 );
200 }
201 }
202 }
203 return iter;
204}
205
206
207Foam::label Foam::snappyRefineDriver::smallFeatureRefine
208(
209 const refinementParameters& refineParams,
210 const label maxIter
211)
212{
213 if (dryRun_)
214 {
215 return 0;
216 }
217
218 if (refineParams.minRefineCells() == -1)
219 {
220 // Special setting to be able to restart shm on meshes with inconsistent
221 // cellLevel/pointLevel
222 return 0;
223 }
224
225 addProfiling(feature, "snappyHexMesh::refine::smallFeature");
226 const fvMesh& mesh = meshRefiner_.mesh();
227
228 label iter = 0;
229
230 // See if any surface has an extendedGapLevel
231 const labelList surfaceMaxLevel(meshRefiner_.surfaces().maxGapLevel());
232 const labelList shellMaxLevel(meshRefiner_.shells().maxGapLevel());
233 const labelList curvMaxLevel(meshRefiner_.surfaces().maxCurvatureLevel());
234
235 if
236 (
237 max(surfaceMaxLevel) == 0
238 && max(shellMaxLevel) == 0
239 && max(curvMaxLevel) == 0
240 )
241 {
242 return iter;
243 }
244
245 for (; iter < maxIter; iter++)
246 {
247 Info<< nl
248 << "Small surface feature refinement iteration " << iter << nl
249 << "--------------------------------------------" << nl
250 << endl;
251
252
253 // Determine cells to refine
254 // ~~~~~~~~~~~~~~~~~~~~~~~~~
255
256 labelList candidateCells
257 (
258 meshRefiner_.refineCandidates
259 (
260 refineParams.locationsInMesh(),
261 refineParams.curvature(),
262 refineParams.planarAngle(),
263
264 false, // featureRefinement
265 false, // featureDistanceRefinement
266 false, // internalRefinement
267 false, // surfaceRefinement
268 false, // curvatureRefinement
269 true, // smallFeatureRefinement
270 false, // gapRefinement
271 false, // bigGapRefinement
272 false, // spreadGapSize
273 refineParams.maxGlobalCells(),
274 refineParams.maxLocalCells()
275 )
276 );
277
278 labelList cellsToRefine
279 (
280 meshRefiner_.meshCutter().consistentRefinement
281 (
282 candidateCells,
283 true
284 )
285 );
286 Info<< "Determined cells to refine in = "
287 << mesh.time().cpuTimeIncrement() << " s" << endl;
288
289
290 const label nCellsToRefine =
291 returnReduce(cellsToRefine.size(), sumOp<label>());
292
293 Info<< "Selected for refinement : " << nCellsToRefine
294 << " cells (out of " << mesh.globalData().nTotalCells()
295 << ')' << endl;
296
297 // Stop when no cells to refine or have done minimum necessary
298 // iterations and not enough cells to refine.
299 if (nCellsToRefine == 0)
300 {
301 Info<< "Stopping refining since too few cells selected."
302 << nl << endl;
303 break;
304 }
305
306
307 if (debug)
308 {
309 const_cast<Time&>(mesh.time())++;
310 }
311
312
313 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
314 {
315 meshRefiner_.balanceAndRefine
316 (
317 "small feature refinement iteration " + name(iter),
318 decomposer_,
319 distributor_,
320 cellsToRefine,
321 refineParams.maxLoadUnbalance(),
322 refineParams.maxCellUnbalance()
323 );
324 }
325 else
326 {
327 meshRefiner_.refineAndBalance
328 (
329 "small feature refinement iteration " + name(iter),
330 decomposer_,
331 distributor_,
332 cellsToRefine,
333 refineParams.maxLoadUnbalance(),
334 refineParams.maxCellUnbalance()
335 );
336 }
337 }
338 return iter;
339}
340
341
342Foam::label Foam::snappyRefineDriver::surfaceOnlyRefine
343(
344 const refinementParameters& refineParams,
345 const label maxIter,
346 const label leakBlockageIter
347)
348{
349 if (dryRun_)
350 {
351 return 0;
352 }
353
354 if (refineParams.minRefineCells() == -1)
355 {
356 // Special setting to be able to restart shm on meshes with inconsistent
357 // cellLevel/pointLevel
358 return 0;
359 }
360
361 addProfiling(surface, "snappyHexMesh::refine::surface");
362 const fvMesh& mesh = meshRefiner_.mesh();
363 const refinementSurfaces& surfaces = meshRefiner_.surfaces();
364
365 // Determine the maximum refinement level over all surfaces. This
366 // determines the minimum number of surface refinement iterations.
367 label overallMaxLevel = max(meshRefiner_.surfaces().maxLevel());
368
369 label iter;
370 for (iter = 0; iter < maxIter; iter++)
371 {
372 Info<< nl
373 << "Surface refinement iteration " << iter << nl
374 << "------------------------------" << nl
375 << endl;
376
377
378 // Do optional leak closing (by removing cells)
379 if (iter >= leakBlockageIter)
380 {
381 // Block off intersections with unzoned surfaces with specified
382 // leakLevel < iter
383 const labelList unnamedSurfaces
384 (
386 (
387 surfaces.surfZones()
388 )
389 );
390
391 DynamicList<label> selectedSurfaces(unnamedSurfaces.size());
392 for (const label surfi : unnamedSurfaces)
393 {
394 const label regioni = surfaces.globalRegion(surfi, 0);
395
396 // Take shortcut: assume all cells on surface are refined to
397 // its refinement level at iteration iter. So just use the
398 // iteration to see if the surface is a candidate.
399 if (iter > surfaces.leakLevel()[regioni])
400 {
401 selectedSurfaces.append(surfi);
402 }
403 }
404
405 if
406 (
407 selectedSurfaces.size()
408 && refineParams.locationsOutsideMesh().size()
409 )
410 {
411 meshRefiner_.blockLeakFaces
412 (
413 globalToMasterPatch_,
414 globalToSlavePatch_,
415 refineParams.locationsInMesh(),
416 refineParams.zonesInMesh(),
417 refineParams.locationsOutsideMesh(),
418 selectedSurfaces
419 );
420 }
421 }
422
423
424 // Determine cells to refine
425 // ~~~~~~~~~~~~~~~~~~~~~~~~~
426 // Only look at surface intersections (minLevel and surface curvature),
427 // do not do internal refinement (refinementShells)
428
429 labelList candidateCells
430 (
431 meshRefiner_.refineCandidates
432 (
433 refineParams.locationsInMesh(),
434 refineParams.curvature(),
435 refineParams.planarAngle(),
436
437 false, // featureRefinement
438 false, // featureDistanceRefinement
439 false, // internalRefinement
440 true, // surfaceRefinement
441 true, // curvatureRefinement
442 false, // smallFeatureRefinement
443 false, // gapRefinement
444 false, // bigGapRefinement
445 false, // spreadGapSize
446 refineParams.maxGlobalCells(),
447 refineParams.maxLocalCells()
448 )
449 );
450 labelList cellsToRefine
451 (
452 meshRefiner_.meshCutter().consistentRefinement
453 (
454 candidateCells,
455 true
456 )
457 );
458 Info<< "Determined cells to refine in = "
459 << mesh.time().cpuTimeIncrement() << " s" << endl;
460
461
462 const label nCellsToRefine =
463 returnReduce(cellsToRefine.size(), sumOp<label>());
464
465 Info<< "Selected for refinement : " << nCellsToRefine
466 << " cells (out of " << mesh.globalData().nTotalCells()
467 << ')' << endl;
468
469 // Stop when no cells to refine or have done minimum necessary
470 // iterations and not enough cells to refine.
471 if
472 (
473 nCellsToRefine == 0
474 || (
475 iter >= overallMaxLevel
476 && nCellsToRefine <= refineParams.minRefineCells()
477 )
478 )
479 {
480 Info<< "Stopping refining since too few cells selected."
481 << nl << endl;
482 break;
483 }
484
485
486 if (debug)
487 {
488 const_cast<Time&>(mesh.time())++;
489 }
490
491
492 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
493 {
494 meshRefiner_.balanceAndRefine
495 (
496 "surface refinement iteration " + name(iter),
497 decomposer_,
498 distributor_,
499 cellsToRefine,
500 refineParams.maxLoadUnbalance(),
501 refineParams.maxCellUnbalance()
502 );
503 }
504 else
505 {
506 meshRefiner_.refineAndBalance
507 (
508 "surface refinement iteration " + name(iter),
509 decomposer_,
510 distributor_,
511 cellsToRefine,
512 refineParams.maxLoadUnbalance(),
513 refineParams.maxCellUnbalance()
514 );
515 }
516 }
517 return iter;
518}
519
520
521Foam::label Foam::snappyRefineDriver::gapOnlyRefine
522(
523 const refinementParameters& refineParams,
524 const label maxIter
525)
526{
527 if (dryRun_)
528 {
529 return 0;
530 }
531
532 if (refineParams.minRefineCells() == -1)
533 {
534 // Special setting to be able to restart shm on meshes with inconsistent
535 // cellLevel/pointLevel
536 return 0;
537 }
538
539 const fvMesh& mesh = meshRefiner_.mesh();
540
541 // Determine the maximum refinement level over all surfaces. This
542 // determines the minimum number of surface refinement iterations.
543
544 label maxIncrement = 0;
545 const labelList& maxLevel = meshRefiner_.surfaces().maxLevel();
546 const labelList& gapLevel = meshRefiner_.surfaces().gapLevel();
547
548 forAll(maxLevel, i)
549 {
550 maxIncrement = max(maxIncrement, gapLevel[i]-maxLevel[i]);
551 }
552
553 label iter = 0;
554
555 if (maxIncrement == 0)
556 {
557 return iter;
558 }
559
560 for (iter = 0; iter < maxIter; iter++)
561 {
562 Info<< nl
563 << "Gap refinement iteration " << iter << nl
564 << "--------------------------" << nl
565 << endl;
566
567
568 // Determine cells to refine
569 // ~~~~~~~~~~~~~~~~~~~~~~~~~
570 // Only look at surface intersections (minLevel and surface curvature),
571 // do not do internal refinement (refinementShells)
572
573 labelList candidateCells
574 (
575 meshRefiner_.refineCandidates
576 (
577 refineParams.locationsInMesh(),
578 refineParams.curvature(),
579 refineParams.planarAngle(),
580
581 false, // featureRefinement
582 false, // featureDistanceRefinement
583 false, // internalRefinement
584 false, // surfaceRefinement
585 false, // curvatureRefinement
586 false, // smallFeatureRefinement
587 true, // gapRefinement
588 false, // bigGapRefinement
589 false, // spreadGapSize
590 refineParams.maxGlobalCells(),
591 refineParams.maxLocalCells()
592 )
593 );
594
596 {
597 Pout<< "Writing current mesh to time "
598 << meshRefiner_.timeName() << endl;
599 meshRefiner_.write
600 (
603 (
606 ),
607 mesh.time().path()/meshRefiner_.timeName()
608 );
609 Pout<< "Dumped mesh in = "
610 << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
611
612
613 Pout<< "Dumping " << candidateCells.size()
614 << " cells to cellSet candidateCellsFromGap." << endl;
615 cellSet c(mesh, "candidateCellsFromGap", candidateCells);
616 c.instance() = meshRefiner_.timeName();
617 c.write();
618 }
619
620 // Grow by one layer to make sure we're covering the gap
621 {
622 boolList isCandidateCell(mesh.nCells(), false);
623 forAll(candidateCells, i)
624 {
625 isCandidateCell[candidateCells[i]] = true;
626 }
627
628 for (label i=0; i<1; i++)
629 {
630 boolList newIsCandidateCell(isCandidateCell);
631
632 // Internal faces
633 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
634 {
635 label own = mesh.faceOwner()[facei];
636 label nei = mesh.faceNeighbour()[facei];
637
638 if (isCandidateCell[own] != isCandidateCell[nei])
639 {
640 newIsCandidateCell[own] = true;
641 newIsCandidateCell[nei] = true;
642 }
643 }
644
645 // Get coupled boundary condition values
646 boolList neiIsCandidateCell;
648 (
649 mesh,
650 isCandidateCell,
651 neiIsCandidateCell
652 );
653
654 // Boundary faces
655 for
656 (
657 label facei = mesh.nInternalFaces();
658 facei < mesh.nFaces();
659 facei++
660 )
661 {
662 label own = mesh.faceOwner()[facei];
663 label bFacei = facei-mesh.nInternalFaces();
664
665 if (isCandidateCell[own] != neiIsCandidateCell[bFacei])
666 {
667 newIsCandidateCell[own] = true;
668 }
669 }
670
671 isCandidateCell.transfer(newIsCandidateCell);
672 }
673
674 label n = 0;
675 forAll(isCandidateCell, celli)
676 {
677 if (isCandidateCell[celli])
678 {
679 n++;
680 }
681 }
682 candidateCells.setSize(n);
683 n = 0;
684 forAll(isCandidateCell, celli)
685 {
686 if (isCandidateCell[celli])
687 {
688 candidateCells[n++] = celli;
689 }
690 }
691 }
692
693
695 {
696 Pout<< "Dumping " << candidateCells.size()
697 << " cells to cellSet candidateCellsFromGapPlusBuffer." << endl;
698 cellSet c(mesh, "candidateCellsFromGapPlusBuffer", candidateCells);
699 c.instance() = meshRefiner_.timeName();
700 c.write();
701 }
702
703
704 labelList cellsToRefine
705 (
706 meshRefiner_.meshCutter().consistentRefinement
707 (
708 candidateCells,
709 true
710 )
711 );
712 Info<< "Determined cells to refine in = "
713 << mesh.time().cpuTimeIncrement() << " s" << endl;
714
715
716 const label nCellsToRefine =
717 returnReduce(cellsToRefine.size(), sumOp<label>());
718
719 Info<< "Selected for refinement : " << nCellsToRefine
720 << " cells (out of " << mesh.globalData().nTotalCells()
721 << ')' << endl;
722
723 // Stop when no cells to refine or have done minimum necessary
724 // iterations and not enough cells to refine.
725 if
726 (
727 nCellsToRefine == 0
728 || (
729 iter >= maxIncrement
730 && nCellsToRefine <= refineParams.minRefineCells()
731 )
732 )
733 {
734 Info<< "Stopping refining since too few cells selected."
735 << nl << endl;
736 break;
737 }
738
739
740 if (debug)
741 {
742 const_cast<Time&>(mesh.time())++;
743 }
744
745
746 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
747 {
748 meshRefiner_.balanceAndRefine
749 (
750 "gap refinement iteration " + name(iter),
751 decomposer_,
752 distributor_,
753 cellsToRefine,
754 refineParams.maxLoadUnbalance(),
755 refineParams.maxCellUnbalance()
756 );
757 }
758 else
759 {
760 meshRefiner_.refineAndBalance
761 (
762 "gap refinement iteration " + name(iter),
763 decomposer_,
764 distributor_,
765 cellsToRefine,
766 refineParams.maxLoadUnbalance(),
767 refineParams.maxCellUnbalance()
768 );
769 }
770 }
771 return iter;
772}
773
774
775Foam::label Foam::snappyRefineDriver::surfaceProximityBlock
776(
777 const refinementParameters& refineParams,
778 const label maxIter
779)
780{
781 if (refineParams.minRefineCells() == -1)
782 {
783 // Special setting to be able to restart shm on meshes with inconsistent
784 // cellLevel/pointLevel
785 return 0;
786 }
787
788 fvMesh& mesh = meshRefiner_.mesh();
789
790 if (min(meshRefiner_.surfaces().blockLevel()) == labelMax)
791 {
792 return 0;
793 }
794
795 label iter = 0;
796
797 for (iter = 0; iter < maxIter; iter++)
798 {
799 Info<< nl
800 << "Gap blocking iteration " << iter << nl
801 << "------------------------" << nl
802 << endl;
803
804
805 // Determine cells to block
806 // ~~~~~~~~~~~~~~~~~~~~~~~~
807
808 meshRefiner_.removeGapCells
809 (
810 refineParams.planarAngle(),
811 meshRefiner_.surfaces().blockLevel(),
812 globalToMasterPatch_,
813 refineParams.nFilterIter()
814 );
815
816 if (debug)
817 {
818 const_cast<Time&>(mesh.time())++;
819 Pout<< "Writing gap blocking iteration "
820 << iter << " mesh to time " << meshRefiner_.timeName()
821 << endl;
822 meshRefiner_.write
823 (
826 (
829 ),
830 mesh.time().path()/meshRefiner_.timeName()
831 );
832 }
833 }
834 return iter;
835}
836
837
838Foam::label Foam::snappyRefineDriver::bigGapOnlyRefine
839(
840 const refinementParameters& refineParams,
841 const bool spreadGapSize,
842 const label maxIter
843)
844{
845 if (refineParams.minRefineCells() == -1)
846 {
847 // Special setting to be able to restart shm on meshes with inconsistent
848 // cellLevel/pointLevel
849 return 0;
850 }
851
852 if (dryRun_)
853 {
854 return 0;
855 }
856
857 const fvMesh& mesh = meshRefiner_.mesh();
858
859 label iter = 0;
860
861 // See if any surface has an extendedGapLevel
862 labelList surfaceMaxLevel(meshRefiner_.surfaces().maxGapLevel());
863 labelList shellMaxLevel(meshRefiner_.shells().maxGapLevel());
864
865 label overallMaxLevel(max(max(surfaceMaxLevel), max(shellMaxLevel)));
866
867 if (overallMaxLevel == 0)
868 {
869 return iter;
870 }
871
872
873 for (; iter < maxIter; iter++)
874 {
875 Info<< nl
876 << "Big gap refinement iteration " << iter << nl
877 << "------------------------------" << nl
878 << endl;
879
880
881 // Determine cells to refine
882 // ~~~~~~~~~~~~~~~~~~~~~~~~~
883
884 labelList candidateCells
885 (
886 meshRefiner_.refineCandidates
887 (
888 refineParams.locationsInMesh(),
889 refineParams.curvature(),
890 refineParams.planarAngle(),
891
892 false, // featureRefinement
893 false, // featureDistanceRefinement
894 false, // internalRefinement
895 false, // surfaceRefinement
896 false, // curvatureRefinement
897 false, // smallFeatureRefinement
898 false, // gapRefinement
899 true, // bigGapRefinement
900 spreadGapSize, // spreadGapSize
901 refineParams.maxGlobalCells(),
902 refineParams.maxLocalCells()
903 )
904 );
905
906
908 {
909 Pout<< "Writing current mesh to time "
910 << meshRefiner_.timeName() << endl;
911 meshRefiner_.write
912 (
915 (
918 ),
919 mesh.time().path()/meshRefiner_.timeName()
920 );
921 Pout<< "Dumped mesh in = "
922 << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
923
924 Pout<< "Dumping " << candidateCells.size()
925 << " cells to cellSet candidateCellsFromBigGap." << endl;
926 cellSet c(mesh, "candidateCellsFromBigGap", candidateCells);
927 c.instance() = meshRefiner_.timeName();
928 c.write();
929 }
930
931 labelList cellsToRefine
932 (
933 meshRefiner_.meshCutter().consistentRefinement
934 (
935 candidateCells,
936 true
937 )
938 );
939 Info<< "Determined cells to refine in = "
940 << mesh.time().cpuTimeIncrement() << " s" << endl;
941
942
943 const label nCellsToRefine =
944 returnReduce(cellsToRefine.size(), sumOp<label>());
945
946 Info<< "Selected for refinement : " << nCellsToRefine
947 << " cells (out of " << mesh.globalData().nTotalCells()
948 << ')' << endl;
949
950 // Stop when no cells to refine or have done minimum necessary
951 // iterations and not enough cells to refine.
952 if
953 (
954 nCellsToRefine == 0
955 || (
956 iter >= overallMaxLevel
957 && nCellsToRefine <= refineParams.minRefineCells()
958 )
959 )
960 {
961 Info<< "Stopping refining since too few cells selected."
962 << nl << endl;
963 break;
964 }
965
966
967 if (debug)
968 {
969 const_cast<Time&>(mesh.time())++;
970 }
971
972
973 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
974 {
975 meshRefiner_.balanceAndRefine
976 (
977 "big gap refinement iteration " + name(iter),
978 decomposer_,
979 distributor_,
980 cellsToRefine,
981 refineParams.maxLoadUnbalance(),
982 refineParams.maxCellUnbalance()
983 );
984 }
985 else
986 {
987 meshRefiner_.refineAndBalance
988 (
989 "big gap refinement iteration " + name(iter),
990 decomposer_,
991 distributor_,
992 cellsToRefine,
993 refineParams.maxLoadUnbalance(),
994 refineParams.maxCellUnbalance()
995 );
996 }
997 }
998 return iter;
999}
1000
1001
1002Foam::label Foam::snappyRefineDriver::danglingCellRefine
1003(
1004 const refinementParameters& refineParams,
1005 const label nFaces,
1006 const label maxIter
1007)
1008{
1009 if (refineParams.minRefineCells() == -1)
1010 {
1011 // Special setting to be able to restart shm on meshes with inconsistent
1012 // cellLevel/pointLevel
1013 return 0;
1014 }
1015
1016 if (dryRun_)
1017 {
1018 return 0;
1019 }
1020
1021 addProfiling(dangling, "snappyHexMesh::refine::danglingCell");
1022 const fvMesh& mesh = meshRefiner_.mesh();
1023
1024 label iter;
1025 for (iter = 0; iter < maxIter; iter++)
1026 {
1027 Info<< nl
1028 << "Dangling coarse cells refinement iteration " << iter << nl
1029 << "--------------------------------------------" << nl
1030 << endl;
1031
1032
1033 // Determine cells to refine
1034 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1035
1036 const cellList& cells = mesh.cells();
1038
1039 labelList candidateCells;
1040 {
1041 cellSet candidateCellSet(mesh, "candidateCells", cells.size()/1000);
1042
1043 forAll(cells, celli)
1044 {
1045 label nIntFaces = 0;
1046 for (const label meshFacei : cells[celli])
1047 {
1048 const label patchi = pbm.patchID(meshFacei);
1049
1050 if (patchi < 0 || pbm[patchi].coupled())
1051 {
1052 ++nIntFaces;
1053 }
1054 }
1055
1056 if (nIntFaces == nFaces)
1057 {
1058 candidateCellSet.insert(celli);
1059 }
1060 }
1061
1063 {
1064 Pout<< "Dumping " << candidateCellSet.size()
1065 << " cells to cellSet candidateCellSet." << endl;
1066 candidateCellSet.instance() = meshRefiner_.timeName();
1067 candidateCellSet.write();
1068 }
1069 candidateCells = candidateCellSet.toc();
1070 }
1071
1072
1073
1074 labelList cellsToRefine
1075 (
1076 meshRefiner_.meshCutter().consistentRefinement
1077 (
1078 candidateCells,
1079 true
1080 )
1081 );
1082 Info<< "Determined cells to refine in = "
1083 << mesh.time().cpuTimeIncrement() << " s" << endl;
1084
1085
1086 const label nCellsToRefine =
1087 returnReduce(cellsToRefine.size(), sumOp<label>());
1088
1089 Info<< "Selected for refinement : " << nCellsToRefine
1090 << " cells (out of " << mesh.globalData().nTotalCells()
1091 << ')' << endl;
1092
1093 // Stop when no cells to refine. After a few iterations check if too
1094 // few cells
1095 if
1096 (
1097 nCellsToRefine == 0
1098 || (
1099 iter >= 1
1100 && nCellsToRefine <= refineParams.minRefineCells()
1101 )
1102 )
1103 {
1104 Info<< "Stopping refining since too few cells selected."
1105 << nl << endl;
1106
1107 if (refineParams.balanceAtEnd())
1108 {
1109 Info<< "Final mesh balancing" << endl;
1110
1111 meshRefiner_.balance
1112 (
1113 "",
1114 decomposer_,
1115 distributor_,
1116 cellsToRefine,
1117 0,
1118 -1
1119 );
1120 }
1121
1122 break;
1123 }
1124
1125
1126 if (debug)
1127 {
1128 const_cast<Time&>(mesh.time())++;
1129 }
1130
1131
1132 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
1133 {
1134 meshRefiner_.balanceAndRefine
1135 (
1136 "coarse cell refinement iteration " + name(iter),
1137 decomposer_,
1138 distributor_,
1139 cellsToRefine,
1140 refineParams.maxLoadUnbalance(),
1141 refineParams.maxCellUnbalance()
1142 );
1143 }
1144 else
1145 {
1146 meshRefiner_.refineAndBalance
1147 (
1148 "coarse cell refinement iteration " + name(iter),
1149 decomposer_,
1150 distributor_,
1151 cellsToRefine,
1152 refineParams.maxLoadUnbalance(),
1153 refineParams.maxCellUnbalance()
1154 );
1155 }
1156 }
1157 return iter;
1158}
1159
1160
1161// Detect cells with opposing intersected faces of differing refinement
1162// level and refine them.
1163Foam::label Foam::snappyRefineDriver::refinementInterfaceRefine
1164(
1165 const refinementParameters& refineParams,
1166 const label maxIter
1167)
1168{
1169 if (refineParams.minRefineCells() == -1)
1170 {
1171 // Special setting to be able to restart shm on meshes with inconsistent
1172 // cellLevel/pointLevel
1173 return 0;
1174 }
1175
1176 if (dryRun_)
1177 {
1178 return 0;
1179 }
1180
1181 addProfiling(interface, "snappyHexMesh::refine::transition");
1182 const fvMesh& mesh = meshRefiner_.mesh();
1183
1184 label iter = 0;
1185
1186 if (refineParams.interfaceRefine())
1187 {
1188 for (;iter < maxIter; iter++)
1189 {
1190 Info<< nl
1191 << "Refinement transition refinement iteration " << iter << nl
1192 << "--------------------------------------------" << nl
1193 << endl;
1194
1195 const labelList& surfaceIndex = meshRefiner_.surfaceIndex();
1196 const hexRef8& cutter = meshRefiner_.meshCutter();
1197 const vectorField& fA = mesh.faceAreas();
1198 const labelList& faceOwner = mesh.faceOwner();
1199
1200
1201 // Determine cells to refine
1202 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1203
1204 const cellList& cells = mesh.cells();
1205
1206 labelList candidateCells;
1207 {
1208 // Pass1: pick up cells with differing face level
1209
1210 cellSet transitionCells
1211 (
1212 mesh,
1213 "transitionCells",
1214 cells.size()/100
1215 );
1216
1217 forAll(cells, celli)
1218 {
1219 const cell& cFaces = cells[celli];
1220 label cLevel = cutter.cellLevel()[celli];
1221
1222 forAll(cFaces, cFacei)
1223 {
1224 label facei = cFaces[cFacei];
1225
1226 if (surfaceIndex[facei] != -1)
1227 {
1228 label fLevel = cutter.faceLevel(facei);
1229 if (fLevel != cLevel)
1230 {
1231 transitionCells.insert(celli);
1232 }
1233 }
1234 }
1235 }
1236
1237
1238 cellSet candidateCellSet
1239 (
1240 mesh,
1241 "candidateCells",
1242 cells.size()/1000
1243 );
1244
1245 // Pass2: check for oppositeness
1246
1247 //for (const label celli : transitionCells)
1248 //{
1249 // const cell& cFaces = cells[celli];
1250 // const point& cc = cellCentres[celli];
1251 // const scalar rCVol = pow(cellVolumes[celli], -5.0/3.0);
1252 //
1253 // // Determine principal axes of cell
1254 // symmTensor R(Zero);
1255 //
1256 // forAll(cFaces, i)
1257 // {
1258 // label facei = cFaces[i];
1259 //
1260 // const point& fc = faceCentres[facei];
1261 //
1262 // // Calculate face-pyramid volume
1263 // scalar pyrVol = 1.0/3.0 * fA[facei] & (fc-cc);
1264 //
1265 // if (faceOwner[facei] != celli)
1266 // {
1267 // pyrVol = -pyrVol;
1268 // }
1269 //
1270 // // Calculate face-pyramid centre
1271 // vector pc = (3.0/4.0)*fc + (1.0/4.0)*cc;
1272 //
1273 // R += pyrVol*sqr(pc-cc)*rCVol;
1274 // }
1275 //
1276 // //- MEJ: Problem: truncation errors cause complex evs
1277 // vector lambdas(eigenValues(R));
1278 // const tensor axes(eigenVectors(R, lambdas));
1279 //
1280 //
1281 // // Check if this cell has
1282 // // - opposing sides intersected
1283 // // - which are of different refinement level
1284 // // - plus the inbetween face
1285 //
1286 // labelVector plusFaceLevel(labelVector(-1, -1, -1));
1287 // labelVector minFaceLevel(labelVector(-1, -1, -1));
1288 //
1289 // forAll(cFaces, cFacei)
1290 // {
1291 // label facei = cFaces[cFacei];
1292 //
1293 // if (surfaceIndex[facei] != -1)
1294 // {
1295 // label fLevel = cutter.faceLevel(facei);
1296 //
1297 // // Get outwards pointing normal
1298 // vector n = fA[facei]/mag(fA[facei]);
1299 // if (faceOwner[facei] != celli)
1300 // {
1301 // n = -n;
1302 // }
1303 //
1304 // // What is major direction and sign
1305 // direction cmpt = vector::X;
1306 // scalar maxComp = (n&axes.x());
1307 //
1308 // scalar yComp = (n&axes.y());
1309 // scalar zComp = (n&axes.z());
1310 //
1311 // if (mag(yComp) > mag(maxComp))
1312 // {
1313 // maxComp = yComp;
1314 // cmpt = vector::Y;
1315 // }
1316 //
1317 // if (mag(zComp) > mag(maxComp))
1318 // {
1319 // maxComp = zComp;
1320 // cmpt = vector::Z;
1321 // }
1322 //
1323 // if (maxComp > 0)
1324 // {
1325 // plusFaceLevel[cmpt] = max
1326 // (
1327 // plusFaceLevel[cmpt],
1328 // fLevel
1329 // );
1330 // }
1331 // else
1332 // {
1333 // minFaceLevel[cmpt] = max
1334 // (
1335 // minFaceLevel[cmpt],
1336 // fLevel
1337 // );
1338 // }
1339 // }
1340 // }
1341 //
1342 // // Check if we picked up any opposite differing level
1343 // for (direction dir = 0; dir < vector::nComponents; dir++)
1344 // {
1345 // if
1346 // (
1347 // plusFaceLevel[dir] != -1
1348 // && minFaceLevel[dir] != -1
1349 // && plusFaceLevel[dir] != minFaceLevel[dir]
1350 // )
1351 // {
1352 // candidateCellSet.insert(celli);
1353 // }
1354 // }
1355 //}
1356
1357 const scalar oppositeCos = Foam::cos(degToRad(135.0));
1358
1359 for (const label celli : transitionCells)
1360 {
1361 const cell& cFaces = cells[celli];
1362 label cLevel = cutter.cellLevel()[celli];
1363
1364 // Detect opposite intersection
1365 bool foundOpposite = false;
1366
1367 forAll(cFaces, cFacei)
1368 {
1369 label facei = cFaces[cFacei];
1370
1371 if
1372 (
1373 surfaceIndex[facei] != -1
1374 && cutter.faceLevel(facei) > cLevel
1375 )
1376 {
1377 // Get outwards pointing normal
1378 vector n = fA[facei]/mag(fA[facei]);
1379 if (faceOwner[facei] != celli)
1380 {
1381 n = -n;
1382 }
1383
1384 // Check for any opposite intersection
1385 forAll(cFaces, cFaceI2)
1386 {
1387 label face2i = cFaces[cFaceI2];
1388
1389 if
1390 (
1391 face2i != facei
1392 && surfaceIndex[face2i] != -1
1393 && cutter.faceLevel(face2i) > cLevel
1394 )
1395 {
1396 // Get outwards pointing normal
1397 vector n2 = fA[face2i]/mag(fA[face2i]);
1398 if (faceOwner[face2i] != celli)
1399 {
1400 n2 = -n2;
1401 }
1402
1403
1404 if ((n&n2) < oppositeCos)
1405 {
1406 foundOpposite = true;
1407 break;
1408 }
1409 }
1410 }
1411
1412 if (foundOpposite)
1413 {
1414 break;
1415 }
1416 }
1417 }
1418
1419
1420 if (foundOpposite)
1421 {
1422 candidateCellSet.insert(celli);
1423 }
1424 }
1425
1427 {
1428 Pout<< "Dumping " << candidateCellSet.size()
1429 << " cells to cellSet candidateCellSet." << endl;
1430 candidateCellSet.instance() = meshRefiner_.timeName();
1431 candidateCellSet.write();
1432 }
1433 candidateCells = candidateCellSet.toc();
1434 }
1435
1436
1437
1438 labelList cellsToRefine
1439 (
1440 meshRefiner_.meshCutter().consistentRefinement
1441 (
1442 candidateCells,
1443 true
1444 )
1445 );
1446 Info<< "Determined cells to refine in = "
1447 << mesh.time().cpuTimeIncrement() << " s" << endl;
1448
1449
1450 const label nCellsToRefine =
1451 returnReduce(cellsToRefine.size(), sumOp<label>());
1452
1453 Info<< "Selected for refinement : " << nCellsToRefine
1454 << " cells (out of " << mesh.globalData().nTotalCells()
1455 << ')' << endl;
1456
1457 // Stop when no cells to refine. After a few iterations check if too
1458 // few cells
1459 if
1460 (
1461 nCellsToRefine == 0
1462 || (
1463 iter >= 1
1464 && nCellsToRefine <= refineParams.minRefineCells()
1465 )
1466 )
1467 {
1468 Info<< "Stopping refining since too few cells selected."
1469 << nl << endl;
1470 break;
1471 }
1472
1473
1474 if (debug)
1475 {
1476 const_cast<Time&>(mesh.time())++;
1477 }
1478
1479
1480 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
1481 {
1482 meshRefiner_.balanceAndRefine
1483 (
1484 "interface cell refinement iteration " + name(iter),
1485 decomposer_,
1486 distributor_,
1487 cellsToRefine,
1488 refineParams.maxLoadUnbalance(),
1489 refineParams.maxCellUnbalance()
1490 );
1491 }
1492 else
1493 {
1494 meshRefiner_.refineAndBalance
1495 (
1496 "interface cell refinement iteration " + name(iter),
1497 decomposer_,
1498 distributor_,
1499 cellsToRefine,
1500 refineParams.maxLoadUnbalance(),
1501 refineParams.maxCellUnbalance()
1502 );
1503 }
1504 }
1505 }
1506 return iter;
1507}
1508
1509bool Foam::snappyRefineDriver::usesHigherLevel
1510(
1511 const labelUList& boundaryPointLevel,
1512 const labelUList& f,
1513 const label cLevel
1514) const
1515{
1516 for (const label pointi : f)
1517 {
1518 if (boundaryPointLevel[pointi] > cLevel)
1519 {
1520 return true;
1521 }
1522 }
1523 return false;
1524}
1525
1526
1527Foam::label Foam::snappyRefineDriver::boundaryRefinementInterfaceRefine
1528(
1529 const refinementParameters& refineParams,
1530 const label maxIter
1531)
1532{
1533 if (refineParams.minRefineCells() == -1)
1534 {
1535 // Special setting to be able to restart shm on meshes with inconsistent
1536 // cellLevel/pointLevel
1537 return 0;
1538 }
1539
1540 if (dryRun_)
1541 {
1542 return 0;
1543 }
1544
1545 addProfiling(interface, "snappyHexMesh::refine::transition");
1546 const fvMesh& mesh = meshRefiner_.mesh();
1547
1548 label iter = 0;
1549
1550 if (refineParams.interfaceRefine())
1551 {
1552 for (;iter < maxIter; iter++)
1553 {
1554 Info<< nl
1555 << "Boundary refinement iteration " << iter << nl
1556 << "-------------------------------" << nl
1557 << endl;
1558
1559 const labelList& surfaceIndex = meshRefiner_.surfaceIndex();
1560 const hexRef8& cutter = meshRefiner_.meshCutter();
1561 const labelList& cellLevel = cutter.cellLevel();
1562 const faceList& faces = mesh.faces();
1563 const cellList& cells = mesh.cells();
1564
1565
1566 // Determine cells to refine
1567 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1568
1569 // Point/face on boundary
1570 bitSet isBoundaryPoint(mesh.nPoints());
1571 bitSet isBoundaryFace(mesh.nFaces());
1572 {
1573 forAll(surfaceIndex, facei)
1574 {
1575 if (surfaceIndex[facei] != -1)
1576 {
1577 isBoundaryFace.set(facei);
1578 isBoundaryPoint.set(faces[facei]);
1579 }
1580 }
1581 const labelList meshPatchIDs(meshRefiner_.meshedPatches());
1582 for (const label patchi : meshPatchIDs)
1583 {
1584 const polyPatch& pp = mesh.boundaryMesh()[patchi];
1585 forAll(pp, i)
1586 {
1587 isBoundaryFace.set(pp.start()+i);
1588 isBoundaryPoint.set(pp[i]);
1589 }
1590 }
1591
1593 (
1594 mesh,
1595 isBoundaryPoint,
1597 0
1598 );
1599 }
1600
1601 // Mark max boundary face level onto boundary points. All points
1602 // not on a boundary face stay 0.
1603 labelList boundaryPointLevel(mesh.nPoints(), 0);
1604 {
1605 forAll(cells, celli)
1606 {
1607 const cell& cFaces = cells[celli];
1608 const label cLevel = cellLevel[celli];
1609
1610 for (const label facei : cFaces)
1611 {
1612 if (isBoundaryFace(facei))
1613 {
1614 const face& f = faces[facei];
1615 for (const label pointi : f)
1616 {
1617 boundaryPointLevel[pointi] =
1618 max
1619 (
1620 boundaryPointLevel[pointi],
1621 cLevel
1622 );
1623 }
1624 }
1625 }
1626 }
1627
1629 (
1630 mesh,
1631 boundaryPointLevel,
1633 label(0)
1634 );
1635 }
1636
1637
1638 // Detect cells with a point but not face on the boundary
1639 labelList candidateCells;
1640 {
1641 const cellList& cells = mesh.cells();
1642
1643 cellSet candidateCellSet
1644 (
1645 mesh,
1646 "candidateCells",
1647 cells.size()/100
1648 );
1649
1650 forAll(cells, celli)
1651 {
1652 const cell& cFaces = cells[celli];
1653 const label cLevel = cellLevel[celli];
1654
1655 bool isBoundaryCell = false;
1656 for (const label facei : cFaces)
1657 {
1658 if (isBoundaryFace(facei))
1659 {
1660 isBoundaryCell = true;
1661 break;
1662 }
1663 }
1664
1665 if (!isBoundaryCell)
1666 {
1667 for (const label facei : cFaces)
1668 {
1669 const face& f = mesh.faces()[facei];
1670 if (usesHigherLevel(boundaryPointLevel, f, cLevel))
1671 {
1672 candidateCellSet.insert(celli);
1673 break;
1674 }
1675 }
1676 }
1677 }
1679 {
1680 Pout<< "Dumping " << candidateCellSet.size()
1681 << " cells to cellSet candidateCellSet." << endl;
1682 candidateCellSet.instance() = meshRefiner_.timeName();
1683 candidateCellSet.write();
1684 }
1685 candidateCells = candidateCellSet.toc();
1686 }
1687
1688 labelList cellsToRefine
1689 (
1690 meshRefiner_.meshCutter().consistentRefinement
1691 (
1692 candidateCells,
1693 true
1694 )
1695 );
1696 Info<< "Determined cells to refine in = "
1697 << mesh.time().cpuTimeIncrement() << " s" << endl;
1698
1699
1700 const label nCellsToRefine =
1701 returnReduce(cellsToRefine.size(), sumOp<label>());
1702
1703 Info<< "Selected for refinement : " << nCellsToRefine
1704 << " cells (out of " << mesh.globalData().nTotalCells()
1705 << ')' << endl;
1706
1707 // Stop when no cells to refine. After a few iterations check if too
1708 // few cells
1709 if
1710 (
1711 nCellsToRefine == 0
1712 // || (
1713 // iter >= 1
1714 // && nCellsToRefine <= refineParams.minRefineCells()
1715 // )
1716 )
1717 {
1718 Info<< "Stopping refining since too few cells selected."
1719 << nl << endl;
1720 break;
1721 }
1722
1723
1724 if (debug)
1725 {
1726 const_cast<Time&>(mesh.time())++;
1727 }
1728
1729
1730 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
1731 {
1732 meshRefiner_.balanceAndRefine
1733 (
1734 "boundary cell refinement iteration " + name(iter),
1735 decomposer_,
1736 distributor_,
1737 cellsToRefine,
1738 refineParams.maxLoadUnbalance(),
1739 refineParams.maxCellUnbalance()
1740 );
1741 }
1742 else
1743 {
1744 meshRefiner_.refineAndBalance
1745 (
1746 "boundary cell refinement iteration " + name(iter),
1747 decomposer_,
1748 distributor_,
1749 cellsToRefine,
1750 refineParams.maxLoadUnbalance(),
1751 refineParams.maxCellUnbalance()
1752 );
1753 }
1754 }
1755 }
1756 return iter;
1757}
1758
1759
1760void Foam::snappyRefineDriver::removeInsideCells
1761(
1762 const refinementParameters& refineParams,
1763 const label nBufferLayers
1764)
1765{
1766 // Skip if no limitRegion and zero bufferLayers
1767 if (meshRefiner_.limitShells().shells().size() == 0 && nBufferLayers == 0)
1768 {
1769 return;
1770 }
1771
1772 if (dryRun_)
1773 {
1774 return;
1775 }
1776
1777 Info<< nl
1778 << "Removing mesh beyond surface intersections" << nl
1779 << "------------------------------------------" << nl
1780 << endl;
1781
1782 const fvMesh& mesh = meshRefiner_.mesh();
1783
1784 if (debug)
1785 {
1786 const_cast<Time&>(mesh.time())++;
1787 }
1788
1789 // Remove any cells inside limitShells with level -1
1790 if (meshRefiner_.limitShells().shells().size())
1791 {
1792 meshRefiner_.removeLimitShells
1793 (
1794 nBufferLayers,
1795 1,
1796 globalToMasterPatch_,
1797 globalToSlavePatch_,
1798 refineParams.locationsInMesh(),
1799 refineParams.zonesInMesh(),
1800 refineParams.locationsOutsideMesh()
1801 );
1802 }
1803
1804
1805 // Fix any additional (e.g. locationsOutsideMesh). Note: probably not
1806 // necessary.
1807 meshRefiner_.splitMesh
1808 (
1809 nBufferLayers, // nBufferLayers
1810 refineParams.nErodeCellZone(),
1811 globalToMasterPatch_,
1812 globalToSlavePatch_,
1813 refineParams.locationsInMesh(),
1814 refineParams.zonesInMesh(),
1815 refineParams.locationsOutsideMesh(),
1816 !refineParams.useLeakClosure(),
1817 setFormatter_,
1818 surfFormatter_
1819 );
1820
1822 {
1823 const_cast<Time&>(mesh.time())++;
1824
1825 Pout<< "Writing subsetted mesh to time "
1826 << meshRefiner_.timeName() << endl;
1827 meshRefiner_.write
1828 (
1831 (
1834 ),
1835 mesh.time().path()/meshRefiner_.timeName()
1836 );
1837 Pout<< "Dumped mesh in = "
1838 << mesh.time().cpuTimeIncrement() << " s\n" << nl << endl;
1839 }
1840}
1841
1842
1843Foam::label Foam::snappyRefineDriver::shellRefine
1844(
1845 const refinementParameters& refineParams,
1846 const label maxIter
1847)
1848{
1849 if (dryRun_)
1850 {
1851 return 0;
1852 }
1853
1854 if (refineParams.minRefineCells() == -1)
1855 {
1856 // Special setting to be able to restart shm on meshes with inconsistent
1857 // cellLevel/pointLevel
1858 return 0;
1859 }
1860
1861 addProfiling(shell, "snappyHexMesh::refine::shell");
1862 const fvMesh& mesh = meshRefiner_.mesh();
1863
1864 // Mark current boundary faces with 0. Have meshRefiner maintain them.
1865 meshRefiner_.userFaceData().setSize(1);
1866
1867 // mark list to remove any refined faces
1868 meshRefiner_.userFaceData()[0].first() = meshRefinement::REMOVE;
1869 meshRefiner_.userFaceData()[0].second() = ListOps::createWithValue<label>
1870 (
1871 mesh.nFaces(),
1872 meshRefiner_.intersectedFaces(),
1873 0, // set value
1874 -1 // default value
1875 );
1876
1877 // Determine the maximum refinement level over all volume refinement
1878 // regions. This determines the minimum number of shell refinement
1879 // iterations.
1880 label overallMaxShellLevel = meshRefiner_.shells().maxLevel();
1881
1882 label iter;
1883 for (iter = 0; iter < maxIter; iter++)
1884 {
1885 Info<< nl
1886 << "Shell refinement iteration " << iter << nl
1887 << "----------------------------" << nl
1888 << endl;
1889
1890 labelList candidateCells
1891 (
1892 meshRefiner_.refineCandidates
1893 (
1894 refineParams.locationsInMesh(),
1895 refineParams.curvature(),
1896 refineParams.planarAngle(),
1897
1898 false, // featureRefinement
1899 true, // featureDistanceRefinement
1900 true, // internalRefinement
1901 false, // surfaceRefinement
1902 false, // curvatureRefinement
1903 false, // smallFeatureRefinement
1904 false, // gapRefinement
1905 false, // bigGapRefinement
1906 false, // spreadGapSize
1907 refineParams.maxGlobalCells(),
1908 refineParams.maxLocalCells()
1909 )
1910 );
1911
1913 {
1914 Pout<< "Dumping " << candidateCells.size()
1915 << " cells to cellSet candidateCellsFromShells." << endl;
1916
1917 cellSet c(mesh, "candidateCellsFromShells", candidateCells);
1918 c.instance() = meshRefiner_.timeName();
1919 c.write();
1920 }
1921
1922 // Problem choosing starting faces for bufferlayers (bFaces)
1923 // - we can't use the current intersected boundary faces
1924 // (intersectedFaces) since this grows indefinitely
1925 // - if we use 0 faces we don't satisfy bufferLayers from the
1926 // surface.
1927 // - possibly we want to have bFaces only the initial set of faces
1928 // and maintain the list while doing the refinement.
1929 labelList bFaces
1930 (
1931 findIndices(meshRefiner_.userFaceData()[0].second(), 0)
1932 );
1933
1934 //Info<< "Collected boundary faces : "
1935 // << returnReduce(bFaces.size(), sumOp<label>()) << endl;
1936
1937 labelList cellsToRefine;
1938
1939 if (refineParams.nBufferLayers() <= 2)
1940 {
1941 cellsToRefine = meshRefiner_.meshCutter().consistentSlowRefinement
1942 (
1943 refineParams.nBufferLayers(),
1944 candidateCells, // cells to refine
1945 bFaces, // faces for nBufferLayers
1946 1, // point difference
1947 meshRefiner_.intersectedPoints() // points to check
1948 );
1949 }
1950 else
1951 {
1952 cellsToRefine = meshRefiner_.meshCutter().consistentSlowRefinement2
1953 (
1954 refineParams.nBufferLayers(),
1955 candidateCells, // cells to refine
1956 bFaces // faces for nBufferLayers
1957 );
1958 }
1959
1960 Info<< "Determined cells to refine in = "
1961 << mesh.time().cpuTimeIncrement() << " s" << endl;
1962
1963
1964 const label nCellsToRefine =
1965 returnReduce(cellsToRefine.size(), sumOp<label>());
1966
1967 Info<< "Selected for internal refinement : " << nCellsToRefine
1968 << " cells (out of " << mesh.globalData().nTotalCells()
1969 << ')' << endl;
1970
1971 // Stop when no cells to refine or have done minimum necessary
1972 // iterations and not enough cells to refine.
1973 if
1974 (
1975 nCellsToRefine == 0
1976 || (
1977 iter >= overallMaxShellLevel
1978 && nCellsToRefine <= refineParams.minRefineCells()
1979 )
1980 )
1981 {
1982 Info<< "Stopping refining since too few cells selected."
1983 << nl << endl;
1984 break;
1985 }
1986
1987
1988 if (debug)
1989 {
1990 const_cast<Time&>(mesh.time())++;
1991 }
1992
1993 if (returnReduceOr(mesh.nCells() >= refineParams.maxLocalCells()))
1994 {
1995 meshRefiner_.balanceAndRefine
1996 (
1997 "shell refinement iteration " + name(iter),
1998 decomposer_,
1999 distributor_,
2000 cellsToRefine,
2001 refineParams.maxLoadUnbalance(),
2002 refineParams.maxCellUnbalance()
2003 );
2004 }
2005 else
2006 {
2007 meshRefiner_.refineAndBalance
2008 (
2009 "shell refinement iteration " + name(iter),
2010 decomposer_,
2011 distributor_,
2012 cellsToRefine,
2013 refineParams.maxLoadUnbalance(),
2014 refineParams.maxCellUnbalance()
2015 );
2016 }
2017 }
2018 meshRefiner_.userFaceData().clear();
2019
2020 return iter;
2021}
2022
2023
2024Foam::label Foam::snappyRefineDriver::directionalShellRefine
2025(
2026 const refinementParameters& refineParams,
2027 const label maxIter
2028)
2029{
2030 if (dryRun_)
2031 {
2032 return 0;
2033 }
2034
2035 addProfiling(shell, "snappyHexMesh::refine::directionalShell");
2036 const fvMesh& mesh = meshRefiner_.mesh();
2037 const shellSurfaces& shells = meshRefiner_.shells();
2038
2039 labelList& cellLevel =
2040 const_cast<labelIOList&>(meshRefiner_.meshCutter().cellLevel());
2041 labelList& pointLevel =
2042 const_cast<labelIOList&>(meshRefiner_.meshCutter().pointLevel());
2043
2044
2045 // Determine the minimum and maximum cell levels that are candidates for
2046 // directional refinement
2047 const labelPairList dirSelect(shells.directionalSelectLevel());
2048 label overallMinLevel = labelMax;
2049 label overallMaxLevel = labelMin;
2050 forAll(dirSelect, shelli)
2051 {
2052 overallMinLevel = min(dirSelect[shelli].first(), overallMinLevel);
2053 overallMaxLevel = max(dirSelect[shelli].second(), overallMaxLevel);
2054 }
2055
2056 if (overallMinLevel > overallMaxLevel)
2057 {
2058 return 0;
2059 }
2060
2061 // Maintain directional refinement levels
2062 List<labelVector> dirCellLevel(cellLevel.size());
2063 forAll(cellLevel, celli)
2064 {
2065 dirCellLevel[celli] = labelVector::uniform(cellLevel[celli]);
2066 }
2067
2068 label iter;
2069 for (iter = 0; iter < maxIter; iter++)
2070 {
2071 Info<< nl
2072 << "Directional shell refinement iteration " << iter << nl
2073 << "----------------------------------------" << nl
2074 << endl;
2075
2076 label nAllRefine = 0;
2077
2078 for (direction dir = 0; dir < vector::nComponents; dir++)
2079 {
2080 // Select the cells that need to be refined in certain direction:
2081 // - cell inside/outside shell
2082 // - original cellLevel (using mapping) mentioned in levelIncrement
2083 // - dirCellLevel not yet up to cellLevel+levelIncrement
2084
2085
2086 // Extract component of directional level
2087 labelList currentLevel(dirCellLevel.size());
2088 forAll(dirCellLevel, celli)
2089 {
2090 currentLevel[celli] = dirCellLevel[celli][dir];
2091 }
2092
2093 labelList candidateCells
2094 (
2095 meshRefiner_.directionalRefineCandidates
2096 (
2097 refineParams.maxGlobalCells(),
2098 refineParams.maxLocalCells(),
2099 currentLevel,
2100 dir
2101 )
2102 );
2103
2104 // Extend to keep 2:1 ratio
2105 labelList cellsToRefine
2106 (
2107 meshRefiner_.meshCutter().consistentRefinement
2108 (
2109 currentLevel,
2110 candidateCells,
2111 true
2112 )
2113 );
2114
2115 Info<< "Determined cells to refine in = "
2116 << mesh.time().cpuTimeIncrement() << " s" << endl;
2117
2118 const label nCellsToRefine =
2119 returnReduce(cellsToRefine.size(), sumOp<label>());
2120
2121 Info<< "Selected for direction " << vector::componentNames[dir]
2122 << " refinement : " << nCellsToRefine
2123 << " cells (out of " << mesh.globalData().nTotalCells()
2124 << ')' << endl;
2125
2126 nAllRefine += nCellsToRefine;
2127
2128 // Stop when no cells to refine or have done minimum necessary
2129 // iterations and not enough cells to refine.
2130 if (nCellsToRefine > 0)
2131 {
2132 if (debug)
2133 {
2134 const_cast<Time&>(mesh.time())++;
2135 }
2136
2137 const bitSet isRefineCell(mesh.nCells(), cellsToRefine);
2138
2140 (
2141 meshRefiner_.directionalRefine
2142 (
2143 "directional refinement iteration " + name(iter),
2144 dir,
2145 cellsToRefine
2146 )
2147 );
2148
2149 Info<< "Refined mesh in = "
2150 << mesh.time().cpuTimeIncrement() << " s" << endl;
2151
2153 (
2154 map().cellMap(),
2155 labelVector(0, 0, 0),
2156 dirCellLevel
2157 );
2158
2159 // Note: edges will have been split. The points might have
2160 // inherited pointLevel from either side of the edge which
2161 // might not be the same for coupled edges so sync
2163 (
2164 mesh,
2165 pointLevel,
2167 labelMin
2168 );
2169
2170 forAll(map().cellMap(), celli)
2171 {
2172 if (isRefineCell[map().cellMap()[celli]])
2173 {
2174 dirCellLevel[celli][dir]++;
2175 }
2176 }
2177
2178 // Do something with the pointLevel. See discussion about the
2179 // cellLevel. Do we keep min/max ?
2180 forAll(map().pointMap(), pointi)
2181 {
2182 label oldPointi = map().pointMap()[pointi];
2183 if (map().reversePointMap()[oldPointi] != pointi)
2184 {
2185 // Is added point (splitting an edge)
2186 pointLevel[pointi]++;
2187 }
2188 }
2189 }
2190 }
2191
2192
2193 if (nAllRefine == 0)
2194 {
2195 Info<< "Stopping refining since no cells selected."
2196 << nl << endl;
2197 break;
2198 }
2199
2200 meshRefiner_.printMeshInfo
2201 (
2202 debug,
2203 "After directional refinement iteration " + name(iter),
2204 true
2205 );
2206
2208 {
2209 Pout<< "Writing directional refinement iteration "
2210 << iter << " mesh to time " << meshRefiner_.timeName() << endl;
2211 meshRefiner_.write
2212 (
2215 (
2218 ),
2219 mesh.time().path()/meshRefiner_.timeName()
2220 );
2221 }
2222 }
2223
2224 // Adjust cellLevel from dirLevel? As max? Or the min?
2225 // For now: use max. The idea is that if there is a wall
2226 // any directional refinement is likely to be aligned with
2227 // the wall (wall layers) so any snapping/layering would probably
2228 // want to use this highest refinement level.
2229
2230 forAll(cellLevel, celli)
2231 {
2232 cellLevel[celli] = cmptMax(dirCellLevel[celli]);
2233 }
2234
2235 return iter;
2236}
2237
2238
2239void Foam::snappyRefineDriver::mergeAndSmoothRatio
2240(
2241 const scalarList& allSeedPointDist,
2242 const label nSmoothExpansion,
2243 List<Tuple2<scalar, scalar>>& keyAndValue
2244)
2245{
2246 // Merge duplicate distance from coupled locations to get unique
2247 // distances to operate on, do on master
2248 SortableList<scalar> unmergedDist(allSeedPointDist);
2249 DynamicList<scalar> mergedDist;
2250
2251 scalar prevDist = GREAT;
2252 forAll(unmergedDist, i)
2253 {
2254 scalar curDist = unmergedDist[i];
2255 scalar difference = mag(curDist - prevDist);
2256 if (difference > meshRefiner_.mergeDistance())
2257 //if (difference > 0.01)
2258 {
2259 mergedDist.append(curDist);
2260 prevDist = curDist;
2261 }
2262 }
2263
2264 // Sort the unique distances
2265 SortableList<scalar> sortedDist(mergedDist);
2266 labelList indexSet = sortedDist.indices();
2267
2268 // Get updated position starting from original (undistorted) mesh
2269 scalarList seedPointsNewLocation = sortedDist;
2270
2271 scalar initResidual = 0.0;
2272 scalar prevIterResidual = GREAT;
2273
2274 for (label iter = 0; iter < nSmoothExpansion; iter++)
2275 {
2276
2277 // Position based edge averaging algorithm operated on
2278 // all seed plane locations in normalized form.
2279 //
2280 // 0 1 2 3 4 5 6 (edge numbers)
2281 // ---x---x---x---x---x---x---
2282 // 0 1 2 3 4 5 (point numbers)
2283 //
2284 // Average of edge 1-3 in terms of position
2285 // = (point3 - point0)/3
2286 // Keeping points 0-1 frozen, new position of point 2
2287 // = position2 + (average of edge 1-3 as above)
2288 for(label i = 2; i<mergedDist.size()-1; i++)
2289 {
2290 scalar oldX00 = sortedDist[i-2];
2291 scalar oldX1 = sortedDist[i+1];
2292 scalar curX0 = seedPointsNewLocation[i-1];
2293 seedPointsNewLocation[i] = curX0 + (oldX1 - oldX00)/3;
2294 }
2295
2296 const scalarField residual(seedPointsNewLocation-sortedDist);
2297 {
2298 scalar res(sumMag(residual));
2299
2300 if (iter == 0)
2301 {
2302 initResidual = res;
2303 }
2304 res /= initResidual;
2305
2306 if (mag(prevIterResidual - res) < SMALL)
2307 {
2308 if (debug)
2309 {
2310 Pout<< "Converged with iteration " << iter
2311 << " initResidual: " << initResidual
2312 << " final residual : " << res << endl;
2313 }
2314 break;
2315 }
2316 else
2317 {
2318 prevIterResidual = res;
2319 }
2320 }
2321
2322 // Update the field for next iteration, avoid moving points
2323 sortedDist = seedPointsNewLocation;
2324
2325 }
2326
2327 keyAndValue.setSize(mergedDist.size());
2328
2329 forAll(mergedDist, i)
2330 {
2331 keyAndValue[i].first() = mergedDist[i];
2332 label index = indexSet[i];
2333 keyAndValue[i].second() = seedPointsNewLocation[index];
2334 }
2335}
2336
2337
2338Foam::label Foam::snappyRefineDriver::directionalSmooth
2339(
2340 const refinementParameters& refineParams
2341)
2342{
2343 addProfiling(split, "snappyHexMesh::refine::smooth");
2344 Info<< nl
2345 << "Directional expansion ratio smoothing" << nl
2346 << "-------------------------------------" << nl
2347 << endl;
2348
2349 fvMesh& baseMesh = meshRefiner_.mesh();
2350 const searchableSurfaces& geometry = meshRefiner_.surfaces().geometry();
2351 const shellSurfaces& shells = meshRefiner_.shells();
2352
2353 label iter = 0;
2354
2355 forAll(shells.nSmoothExpansion(), shellI)
2356 {
2357 if
2358 (
2359 shells.nSmoothExpansion()[shellI] > 0
2360 || shells.nSmoothPosition()[shellI] > 0
2361 )
2362 {
2363 label surfi = shells.shells()[shellI];
2364 const vector& userDirection = shells.smoothDirection()[shellI];
2365
2366
2367 // Extract inside points
2369 {
2370 // Get inside points
2371 List<volumeType> volType;
2372 geometry[surfi].getVolumeType(baseMesh.points(), volType);
2373
2374 label nInside = 0;
2375 forAll(volType, pointi)
2376 {
2377 if (volType[pointi] == volumeType::INSIDE)
2378 {
2379 nInside++;
2380 }
2381 }
2382 pointLabels.setSize(nInside);
2383 nInside = 0;
2384 forAll(volType, pointi)
2385 {
2386 if (volType[pointi] == volumeType::INSIDE)
2387 {
2388 pointLabels[nInside++] = pointi;
2389 }
2390 }
2391
2392 //bitSet isInsidePoint(baseMesh.nPoints());
2393 //forAll(volType, pointi)
2394 //{
2395 // if (volType[pointi] == volumeType::INSIDE)
2396 // {
2397 // isInsidePoint.set(pointi);
2398 // }
2399 //}
2400 //pointLabels = isInsidePoint.used();
2401 }
2402
2403 // Mark all directed edges
2404 bitSet isXEdge(baseMesh.edges().size());
2405 {
2406 const edgeList& edges = baseMesh.edges();
2407 forAll(edges, edgei)
2408 {
2409 const edge& e = edges[edgei];
2410 vector eVec(e.vec(baseMesh.points()));
2411 eVec /= mag(eVec);
2412 if (mag(eVec&userDirection) > 0.9)
2413 {
2414 isXEdge.set(edgei);
2415 }
2416 }
2417 }
2418
2419 // Get the extreme of smoothing region and
2420 // normalize all points within
2421 const scalar totalLength =
2422 geometry[surfi].bounds().span()
2423 & userDirection;
2424 const scalar startPosition =
2425 geometry[surfi].bounds().min()
2426 & userDirection;
2427
2428 scalarField normalizedPosition(pointLabels.size(), Zero);
2430 {
2431 label pointi = pointLabels[i];
2432 normalizedPosition[i] =
2433 (
2434 ((baseMesh.points()[pointi]&userDirection) - startPosition)
2435 / totalLength
2436 );
2437 }
2438
2439 // Sort the normalized position
2440 labelList order(sortedOrder(normalizedPosition));
2441
2442 DynamicList<scalar> seedPointDist;
2443
2444 // Select points from finest refinement (one point-per plane)
2445 scalar prevDist = GREAT;
2446 forAll(order, i)
2447 {
2448 label pointi = order[i];
2449 scalar curDist = normalizedPosition[pointi];
2450 if (mag(curDist - prevDist) > meshRefiner_.mergeDistance())
2451 {
2452 seedPointDist.append(curDist);
2453 prevDist = curDist;
2454 }
2455 }
2456
2457 // Collect data from all processors
2458 scalarList allSeedPointDist;
2459 {
2460 List<scalarList> gatheredDist(Pstream::nProcs());
2461 gatheredDist[Pstream::myProcNo()] = seedPointDist;
2462 Pstream::gatherList(gatheredDist);
2463
2464 // Combine processor lists into one big list.
2465 allSeedPointDist =
2467 (
2468 gatheredDist, accessOp<scalarList>()
2469 );
2470 }
2471
2472 // Pre-set the points not to smooth (after expansion)
2473 bitSet isFrozenPoint(baseMesh.nPoints(), true);
2474
2475 {
2476 scalar minSeed = min(allSeedPointDist);
2477 scalar maxSeed = max(allSeedPointDist);
2478 Pstream::broadcasts(UPstream::worldComm, minSeed, maxSeed);
2479
2480 forAll(normalizedPosition, posI)
2481 {
2482 const scalar pos = normalizedPosition[posI];
2483 if
2484 (
2485 (mag(pos-minSeed) < meshRefiner_.mergeDistance())
2486 || (mag(pos-maxSeed) < meshRefiner_.mergeDistance())
2487 )
2488 {
2489 // Boundary point: freeze
2490 isFrozenPoint.set(pointLabels[posI]);
2491 }
2492 else
2493 {
2494 // Internal to moving region
2495 isFrozenPoint.unset(pointLabels[posI]);
2496 }
2497 }
2498 }
2499
2500 Info<< "Smoothing " << geometry[surfi].name() << ':' << nl
2501 << " Direction : " << userDirection << nl
2502 << " Number of points : "
2504 << " (out of " << baseMesh.globalData().nTotalPoints()
2505 << ")" << nl
2506 << " Smooth expansion iterations : "
2507 << shells.nSmoothExpansion()[shellI] << nl
2508 << " Smooth position iterations : "
2509 << shells.nSmoothPosition()[shellI] << nl
2510 << " Number of planes : "
2511 << allSeedPointDist.size()
2512 << endl;
2513
2514 // Make lookup from original normalized distance to new value
2515 List<Tuple2<scalar, scalar>> keyAndValue(allSeedPointDist.size());
2516
2517 // Filter unique seed distances and iterate for user given steps
2518 // or until convergence. Then get back map from old to new distance
2519 if (Pstream::master())
2520 {
2521 mergeAndSmoothRatio
2522 (
2523 allSeedPointDist,
2524 shells.nSmoothExpansion()[shellI],
2525 keyAndValue
2526 );
2527 }
2528
2529 Pstream::broadcast(keyAndValue);
2530
2531 // Construct an iterpolation table for further queries
2532 // - although normalized values are used for query,
2533 // it might flow out of bounds due to precision, hence clamped
2534 const interpolationTable<scalar> table
2535 (
2536 keyAndValue,
2538 "undefined"
2539 );
2540
2541 // Now move the points directly on the baseMesh.
2542 pointField baseNewPoints(baseMesh.points());
2544 {
2545 label pointi = pointLabels[i];
2546 const point& curPoint = baseMesh.points()[pointi];
2547 scalar curDist = normalizedPosition[i];
2548 scalar newDist = table(curDist);
2549 scalar newPosition = startPosition + newDist*totalLength;
2550 baseNewPoints[pointi] +=
2551 userDirection * (newPosition - (curPoint &userDirection));
2552 }
2553
2554 // Moving base mesh with expansion ratio smoothing
2555 vectorField disp(baseNewPoints-baseMesh.points());
2557 (
2558 baseMesh,
2559 disp,
2562 );
2563 baseMesh.movePoints(baseMesh.points()+disp);
2564
2565 // Reset any moving state
2566 baseMesh.moving(false);
2567
2569 {
2570 const_cast<Time&>(baseMesh.time())++;
2571
2572 Pout<< "Writing directional expansion ratio smoothed"
2573 << " mesh to time " << meshRefiner_.timeName() << endl;
2574
2575 meshRefiner_.write
2576 (
2579 (
2582 ),
2583 baseMesh.time().path()/meshRefiner_.timeName()
2584 );
2585 }
2586
2587 // Now we have moved the points in user specified region. Smooth
2588 // them with neighbour points to avoid skewed cells
2589 // Instead of moving actual mesh, operate on copy
2590 pointField baseMeshPoints(baseMesh.points());
2591 scalar initResidual = 0.0;
2592 scalar prevIterResidual = GREAT;
2593 for (iter = 0; iter < shells.nSmoothPosition()[shellI]; iter++)
2594 {
2595 {
2596 const edgeList& edges = baseMesh.edges();
2597 const labelListList& pointEdges = baseMesh.pointEdges();
2598
2599 pointField unsmoothedPoints(baseMeshPoints);
2600
2601 scalarField sumOther(baseMesh.nPoints(), Zero);
2602 labelList nSumOther(baseMesh.nPoints(), Zero);
2603 labelList nSumXEdges(baseMesh.nPoints(), Zero);
2604 forAll(edges, edgei)
2605 {
2606 const edge& e = edges[edgei];
2607 sumOther[e[0]] +=
2608 (unsmoothedPoints[e[1]]&userDirection);
2609 nSumOther[e[0]]++;
2610 sumOther[e[1]] +=
2611 (unsmoothedPoints[e[0]]&userDirection);
2612 nSumOther[e[1]]++;
2613 if (isXEdge[edgei])
2614 {
2615 nSumXEdges[e[0]]++;
2616 nSumXEdges[e[1]]++;
2617 }
2618 }
2619
2621 (
2622 baseMesh,
2623 nSumXEdges,
2625 label(0)
2626 );
2627
2629 {
2630 label pointi = pointLabels[i];
2631
2632 if (nSumXEdges[pointi] < 2)
2633 {
2634 // Hanging node. Remove the (single!) X edge so it
2635 // will follow points above or below instead
2636 const labelList& pEdges = pointEdges[pointi];
2637 forAll(pEdges, pE)
2638 {
2639 label edgei = pEdges[pE];
2640 if (isXEdge[edgei])
2641 {
2642 const edge& e = edges[edgei];
2643 label otherPt = e.otherVertex(pointi);
2644 nSumOther[pointi]--;
2645 sumOther[pointi] -=
2646 (
2647 unsmoothedPoints[otherPt]
2648 & userDirection
2649 );
2650 }
2651 }
2652 }
2653 }
2654
2656 (
2657 baseMesh,
2658 sumOther,
2660 scalar(0)
2661 );
2663 (
2664 baseMesh,
2665 nSumOther,
2667 label(0)
2668 );
2669
2671 {
2672 label pointi = pointLabels[i];
2673
2674 if ((nSumOther[pointi] >= 2) && !isFrozenPoint[pointi])
2675 {
2676 scalar smoothPos =
2677 0.5
2678 *(
2679 (unsmoothedPoints[pointi]&userDirection)
2680 +sumOther[pointi]/nSumOther[pointi]
2681 );
2682
2683 vector& v = baseNewPoints[pointi];
2684 v += (smoothPos-(v&userDirection))*userDirection;
2685 }
2686 }
2687
2688 const vectorField residual(baseNewPoints - baseMeshPoints);
2689 {
2690 scalar res(gSumMag(residual));
2691
2692 if (iter == 0)
2693 {
2694 initResidual = res;
2695 }
2696 res /= initResidual;
2697
2698 if (mag(prevIterResidual - res) < SMALL)
2699 {
2700 Info<< "Converged smoothing in iteration " << iter
2701 << " initResidual: " << initResidual
2702 << " final residual : " << res << endl;
2703 break;
2704 }
2705 else
2706 {
2707 prevIterResidual = res;
2708 }
2709 }
2710
2711 // Just copy new location instead of moving base mesh
2712 baseMeshPoints = baseNewPoints;
2713 }
2714 }
2715
2716 // Move mesh to new location
2717 vectorField dispSmooth(baseMeshPoints-baseMesh.points());
2719 (
2720 baseMesh,
2721 dispSmooth,
2724 );
2725 baseMesh.movePoints(baseMesh.points()+dispSmooth);
2726
2727 // Reset any moving state
2728 baseMesh.moving(false);
2729
2731 {
2732 const_cast<Time&>(baseMesh.time())++;
2733
2734 Pout<< "Writing positional smoothing iteration "
2735 << iter << " mesh to time " << meshRefiner_.timeName()
2736 << endl;
2737 meshRefiner_.write
2738 (
2741 (
2744 ),
2745 baseMesh.time().path()/meshRefiner_.timeName()
2746 );
2747 }
2748 }
2749 }
2750 return iter;
2751}
2752
2753
2754void Foam::snappyRefineDriver::baffleAndSplitMesh
2755(
2756 const refinementParameters& refineParams,
2757 const snapParameters& snapParams,
2758 const bool handleSnapProblems,
2759 const dictionary& motionDict
2760)
2761{
2762 if (dryRun_)
2763 {
2764 return;
2765 }
2766
2767 addProfiling(split, "snappyHexMesh::refine::splitting");
2768 Info<< nl
2769 << "Splitting mesh at surface intersections" << nl
2770 << "---------------------------------------" << nl
2771 << endl;
2772
2773 const fvMesh& mesh = meshRefiner_.mesh();
2774
2775 if (debug)
2776 {
2777 const_cast<Time&>(mesh.time())++;
2778 }
2779
2780 // Introduce baffles at surface intersections. Note:
2781 // meshRefinement::surfaceIndex() will
2782 // be like boundary face from now on so not coupled anymore.
2783 meshRefiner_.baffleAndSplitMesh
2784 (
2785 handleSnapProblems, // detect&remove potential snap problem
2786
2787 // Snap problem cell detection
2788 snapParams,
2789 refineParams.useTopologicalSnapDetection(),
2790 false, // perpendicular edge connected cells
2791 scalarField(0), // per region perpendicular angle
2792 refineParams.nErodeCellZone(),
2793
2794 motionDict,
2795 const_cast<Time&>(mesh.time()),
2796 globalToMasterPatch_,
2797 globalToSlavePatch_,
2798 refineParams.locationsInMesh(),
2799 refineParams.zonesInMesh(),
2800 refineParams.locationsOutsideMesh(),
2801 !refineParams.useLeakClosure(),
2802 setFormatter_,
2803 surfFormatter_
2804 );
2805
2806
2807 if (!handleSnapProblems) // merge free standing baffles?
2808 {
2809 // By default only merge baffles if on same patch. However the patching
2810 // of the castellated mesh is not very accurate. For backwards
2811 // compatibility disable check for new mode only.
2812 const auto mt = meshRefiner_.meshType();
2813 const bool mergeSameOnly
2814 (
2815 (
2818 )
2819 ? false
2820 : true
2821 );
2822
2823 meshRefiner_.mergeFreeStandingBaffles
2824 (
2825 mergeSameOnly, // whether to merge baffles on different patches
2826 snapParams,
2827 refineParams.useTopologicalSnapDetection(),
2828 false, // perpendicular edge connected cells
2829 scalarField(0), // per region perpendicular angle
2830 refineParams.planarAngle(),
2831 motionDict,
2832 const_cast<Time&>(mesh.time()),
2833 globalToMasterPatch_,
2834 globalToSlavePatch_,
2835 refineParams.locationsInMesh(),
2836 refineParams.locationsOutsideMesh()
2837 );
2838 }
2839}
2840
2841
2842void Foam::snappyRefineDriver::zonify
2843(
2844 const refinementParameters& refineParams,
2845 wordPairHashTable& zonesToFaceZone
2846)
2847{
2848 // Mesh is at its finest. Do zoning
2849 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2850 // This puts all faces with intersection across a zoneable surface
2851 // into that surface's faceZone. All cells inside faceZone get given the
2852 // same cellZone.
2853
2854 const labelList namedSurfaces =
2855 surfaceZonesInfo::getNamedSurfaces(meshRefiner_.surfaces().surfZones());
2856
2857 if
2858 (
2859 namedSurfaces.size()
2860 || refineParams.zonesInMesh().size()
2861 )
2862 {
2863 Info<< nl
2864 << "Introducing zones for interfaces" << nl
2865 << "--------------------------------" << nl
2866 << endl;
2867
2868 const fvMesh& mesh = meshRefiner_.mesh();
2869
2870 if (debug)
2871 {
2872 const_cast<Time&>(mesh.time())++;
2873 }
2874
2875 meshRefiner_.zonify
2876 (
2877 refineParams.allowFreeStandingZoneFaces(),
2878 refineParams.nErodeCellZone(),
2879 refineParams.locationsInMesh(),
2880 refineParams.zonesInMesh(),
2881 refineParams.locationsOutsideMesh(),
2882 !refineParams.useLeakClosure(),
2883 setFormatter_,
2884 surfFormatter_,
2885 zonesToFaceZone
2886 );
2887
2889 {
2890 Pout<< "Writing zoned mesh to time "
2891 << meshRefiner_.timeName() << endl;
2892 meshRefiner_.write
2893 (
2896 (
2899 ),
2900 mesh.time().path()/meshRefiner_.timeName()
2901 );
2902 }
2903
2904 // Check that all faces are synced
2906 }
2907}
2908
2909
2910void Foam::snappyRefineDriver::splitAndMergeBaffles
2911(
2912 const refinementParameters& refineParams,
2913 const snapParameters& snapParams,
2914 const bool handleSnapProblems,
2915 const dictionary& motionDict
2916)
2917{
2918 if (dryRun_)
2919 {
2920 return;
2921 }
2922
2923 Info<< nl
2924 << "Handling cells with snap problems" << nl
2925 << "---------------------------------" << nl
2926 << endl;
2927
2928 const fvMesh& mesh = meshRefiner_.mesh();
2929
2930 // Introduce baffles and split mesh
2931 if (debug)
2932 {
2933 const_cast<Time&>(mesh.time())++;
2934 }
2935
2936 const scalarField& perpAngle = meshRefiner_.surfaces().perpendicularAngle();
2937
2938 meshRefiner_.baffleAndSplitMesh
2939 (
2940 handleSnapProblems,
2941
2942 // Snap problem cell detection
2943 snapParams,
2944 refineParams.useTopologicalSnapDetection(),
2945 handleSnapProblems, // remove perp edge connected cells
2946 perpAngle, // perp angle
2947 refineParams.nErodeCellZone(),
2948
2949 motionDict,
2950 const_cast<Time&>(mesh.time()),
2951 globalToMasterPatch_,
2952 globalToSlavePatch_,
2953 refineParams.locationsInMesh(),
2954 refineParams.zonesInMesh(),
2955 refineParams.locationsOutsideMesh(),
2956 !refineParams.useLeakClosure(),
2957 setFormatter_,
2958 surfFormatter_
2959 );
2960
2961
2962 // By default only merge baffles if on same patch. However the patching
2963 // of the castellated mesh is not very accurate. For backwards
2964 // compatibility disable check for new mode only.
2965 const auto mt = meshRefiner_.meshType();
2966 const bool mergeSameOnly
2967 (
2968 (
2971 )
2972 ? false
2973 : true
2974 );
2975
2976 meshRefiner_.mergeFreeStandingBaffles
2977 (
2978 mergeSameOnly,
2979 snapParams,
2980 refineParams.useTopologicalSnapDetection(),
2981 handleSnapProblems,
2982 perpAngle,
2983 refineParams.planarAngle(),
2984 motionDict,
2985 const_cast<Time&>(mesh.time()),
2986 globalToMasterPatch_,
2987 globalToSlavePatch_,
2988 refineParams.locationsInMesh(),
2989 refineParams.locationsOutsideMesh()
2990 );
2991
2992 if (debug)
2993 {
2994 const_cast<Time&>(mesh.time())++;
2995 }
2996
2997 // Duplicate points on baffles that are on more than one cell
2998 // region. This will help snapping pull them to separate surfaces.
2999 meshRefiner_.dupNonManifoldPoints();
3000
3001
3002 // Merge all baffles that are still remaining after duplicating points.
3004
3005 const label nCouples = returnReduce(couples.size(), sumOp<label>());
3006
3007 Info<< "Detected unsplittable baffles : " << nCouples << endl;
3008
3009 if (nCouples > 0)
3010 {
3011 // Actually merge baffles. Note: not exactly parallellized. Should
3012 // convert baffle faces into processor faces if they resulted
3013 // from them.
3014 meshRefiner_.mergeBaffles(couples, Map<label>(0));
3015
3016 if (debug)
3017 {
3018 // Debug:test all is still synced across proc patches
3019 meshRefiner_.checkData();
3020 }
3021
3022 // Remove any now dangling parts
3023 meshRefiner_.splitMeshRegions
3024 (
3025 globalToMasterPatch_,
3026 globalToSlavePatch_,
3027 refineParams.locationsInMesh(),
3028 refineParams.locationsOutsideMesh(),
3029 true, // exit if any path to outside points
3030 setFormatter_
3031 );
3032
3033 if (debug)
3034 {
3035 // Debug:test all is still synced across proc patches
3036 meshRefiner_.checkData();
3037 }
3038
3039 Info<< "Merged free-standing baffles in = "
3040 << mesh.time().cpuTimeIncrement() << " s." << endl;
3041 }
3042
3044 {
3045 Pout<< "Writing handleProblemCells mesh to time "
3046 << meshRefiner_.timeName() << endl;
3047 meshRefiner_.write
3048 (
3051 (
3054 ),
3055 mesh.time().path()/meshRefiner_.timeName()
3056 );
3057 }
3058}
3059
3060
3061void Foam::snappyRefineDriver::erodeNonManifoldZoneFaces
3062(
3063 const refinementParameters& refineParams
3064)
3065{
3066 if (dryRun_)
3067 {
3068 return;
3069 }
3070
3071 addProfiling(merge, "snappyHexMesh::refine::erode");
3072 Info<< nl
3073 << "Erode non-manifold zone faces" << nl
3074 << "-----------------------------" << nl
3075 << endl;
3076
3077 auto& mesh = meshRefiner_.mesh();
3078 auto& fzs = mesh.faceZones();
3079
3080 // Detect a faceZone face that:
3081 // - sits on the edge of the faceZone (i.e. has an open edge)
3082 // - has a non-manifold edge
3083 //
3084 // Limitations:
3085 // - analyses across all zones - does not support overlapping zones
3086 // - requires zones to be consistent across processor patches
3087
3088 bitSet isZoneFace(mesh.nFaces());
3089 for (const auto& fz : fzs)
3090 {
3091 isZoneFace.set(fz.addressing());
3092 }
3093 // Unmark non-owner side of processor patches so they don't get double
3094 // counted.
3095 const auto& pbm = mesh.boundaryMesh();
3096 for (label patchi = pbm.nNonProcessor(); patchi < pbm.size(); patchi++)
3097 {
3098 if (!refCast<const processorPolyPatch>(pbm[patchi]).owner())
3099 {
3100 isZoneFace.unset(pbm[patchi].range());
3101 }
3102 }
3103
3104
3105 // TBD: replace with mesh.globalData().globalMeshFaceAddr()
3106 const globalIndex globalFaces(mesh.nFaces());
3107
3108 label nChanged = 0;
3109
3110 while (true)
3111 {
3112 const labelList meshFaces(isZoneFace.sortedToc());
3113
3115 (
3117 (
3118 mesh.faces(),
3119 meshFaces
3120 ),
3121 mesh.points()
3122 );
3123
3124 const labelListList edgeGlobalFaces
3125 (
3127 (
3128 mesh,
3129 globalFaces,
3130 pp
3131 )
3132 );
3133
3134 const label nOldChanged = nChanged;
3135 forAll(pp, patchFacei)
3136 {
3137 // Detect at least one open edge and one non-manifold edge
3138 label nOpen = 0;
3139 label nNonManif = 0;
3140 for (const label edgei : pp.faceEdges()[patchFacei])
3141 {
3142 if (edgeGlobalFaces[edgei].size() == 1)
3143 {
3144 nOpen++;
3145 }
3146 else if (edgeGlobalFaces[edgei].size() > 2)
3147 {
3148 nNonManif++;
3149 }
3150 }
3151 if (nNonManif > 0 && nOpen > 0)
3152 {
3153 isZoneFace.unset(meshFaces[patchFacei]);
3154 nChanged++;
3155 }
3156 }
3157
3158 if (returnReduce((nChanged-nOldChanged), sumOp<label>()) == 0)
3159 {
3160 break;
3161 }
3162 }
3163
3164
3165 reduce(nChanged, sumOp<label>());
3166
3167 if (nChanged)
3168 {
3169 // Filter out eroded faceZone faces
3170 fzs.clearAddressing();
3171
3172 for (auto& fz : fzs)
3173 {
3174 DynamicList<label> addressing(fz.size());
3175 DynamicList<bool> flipMap(fz.size());
3176
3177 forAll(fz.addressing(), i)
3178 {
3179 if (isZoneFace[fz.addressing()[i]])
3180 {
3181 addressing.append(fz.addressing()[i]);
3182 flipMap.append(fz.flipMap()[i]);
3183 }
3184 }
3185
3186 //Pout<< "** on faceZone:" << fz.name()
3187 // << " from:" << fz.size()
3188 // << " to:" << addressing.size() << endl;
3189
3190 fz.resetAddressing(addressing, flipMap);
3191 }
3192 }
3194 Info<< "Eroded " << nChanged
3195 << " free-standing zone faces in = "
3196 << mesh.time().cpuTimeIncrement() << " s." << endl;
3197}
3198
3199
3201(
3202 meshRefinement& meshRefiner,
3203 const refinementParameters& refineParams,
3204 const HashTable<Pair<word>>& faceZoneToPatches
3205)
3206{
3207 if (faceZoneToPatches.size())
3208 {
3209 Info<< nl
3210 << "Adding patches for face zones" << nl
3211 << "-----------------------------" << nl
3212 << endl;
3213
3214 Info<< setf(ios_base::left)
3215 << setw(6) << "Patch"
3216 << setw(20) << "Type"
3217 << setw(30) << "Name"
3218 << setw(30) << "FaceZone"
3219 << setw(10) << "FaceType"
3220 << nl
3221 << setw(6) << "-----"
3222 << setw(20) << "----"
3223 << setw(30) << "----"
3224 << setw(30) << "--------"
3225 << setw(10) << "--------"
3226 << endl;
3227
3228 const polyMesh& mesh = meshRefiner.mesh();
3229
3230 // Add patches for added inter-region faceZones
3231 forAllConstIters(faceZoneToPatches, iter)
3232 {
3233 const word& fzName = iter.key();
3234 const Pair<word>& patchNames = iter.val();
3235
3236 // Get any user-defined faceZone data
3238 dictionary patchInfo = refineParams.getZoneInfo(fzName, fzType);
3239
3240 const word& masterName = fzName;
3241 //const word slaveName = fzName + "_slave";
3242 //const word slaveName = czNames.second()+"_to_"+czNames.first();
3243 const word& slaveName = patchNames.second();
3244
3245 label mpi = meshRefiner.addMeshedPatch(masterName, patchInfo);
3246
3247 Info<< setf(ios_base::left)
3248 << setw(6) << mpi
3249 << setw(20) << mesh.boundaryMesh()[mpi].type()
3250 << setw(30) << masterName
3251 << setw(30) << fzName
3253 << nl;
3254
3255
3256 label sli = meshRefiner.addMeshedPatch(slaveName, patchInfo);
3257
3258 Info<< setf(ios_base::left)
3259 << setw(6) << sli
3260 << setw(20) << mesh.boundaryMesh()[sli].type()
3261 << setw(30) << slaveName
3262 << setw(30) << fzName
3264 << nl;
3265
3266 meshRefiner.addFaceZone(fzName, masterName, slaveName, fzType);
3267 }
3268
3269 Info<< endl;
3270 }
3271}
3272
3273
3274void Foam::snappyRefineDriver::mergePatchFaces
3275(
3276 const meshRefinement::FaceMergeType mergeType,
3277 const refinementParameters& refineParams,
3278 const dictionary& motionDict
3279)
3280{
3281 if (dryRun_)
3282 {
3283 return;
3284 }
3285
3286 addProfiling(merge, "snappyHexMesh::refine::merge");
3287 Info<< nl
3288 << "Merge refined boundary faces" << nl
3289 << "----------------------------" << nl
3290 << endl;
3291
3292 const fvMesh& mesh = meshRefiner_.mesh();
3293
3294 if
3295 (
3298 )
3299 {
3300 meshRefiner_.mergePatchFacesUndo
3301 (
3302 Foam::cos(degToRad(45.0)),
3303 Foam::cos(degToRad(45.0)),
3304 meshRefiner_.meshedPatches(),
3305 motionDict,
3306 labelList(mesh.nFaces(), -1),
3307 mergeType
3308 );
3309 }
3310 else
3311 {
3312 // Still merge refined boundary faces if all four are on same patch
3313 meshRefiner_.mergePatchFaces
3314 (
3315 Foam::cos(degToRad(45.0)),
3316 Foam::cos(degToRad(45.0)),
3317 4, // only merge faces split into 4
3318 meshRefiner_.meshedPatches(),
3319 meshRefinement::FaceMergeType::GEOMETRIC // no merge across patches
3320 );
3321 }
3322
3323 if (debug)
3324 {
3325 meshRefiner_.checkData();
3326 }
3327
3328 meshRefiner_.mergeEdgesUndo(Foam::cos(degToRad(45.0)), motionDict);
3329
3330 if (debug)
3331 {
3332 meshRefiner_.checkData();
3333 }
3334}
3335
3336
3337void Foam::snappyRefineDriver::deleteSmallRegions
3338(
3339 const refinementParameters& refineParams
3340)
3341{
3342 const fvMesh& mesh = meshRefiner_.mesh();
3343 const cellZoneMesh& czm = mesh.cellZones();
3344
3345 //const labelList zoneIDs
3346 //(
3347 // meshRefiner_.getZones
3348 // (
3349 // List<surfaceZonesInfo::faceZoneType> fzTypes
3350 // ({
3351 // surfaceZonesInfo::BAFFLE,
3352 // surfaceZonesInfo::BOUNDARY,
3353 // });
3354 // )
3355 //);
3357
3358 // Get faceZone and patch(es) per face (or -1 if face not on faceZone)
3360 labelList ownPatch;
3361 labelList neiPatch;
3362 labelList nB; // local count of faces per faceZone
3363 meshRefiner_.getZoneFaces(zoneIDs, faceZoneID, ownPatch, neiPatch, nB);
3364
3365
3366 // Mark all faces on outside of zones. Note: assumes that faceZones
3367 // are consistent with the outside of cellZones ...
3368
3370 meshRefiner_.selectSeparatedCoupledFaces(isBlockedFace);
3371
3372 forAll(ownPatch, facei)
3373 {
3374 if (ownPatch[facei] != -1)
3375 {
3376 isBlockedFace[facei] = true;
3377 }
3378 }
3379
3380 // Map from cell to zone. Not that background cells are not -1 but
3381 // cellZones.size()
3382 labelList cellToZone(mesh.nCells(), czm.size());
3383 for (const auto& cz : czm)
3384 {
3385 UIndirectList<label>(cellToZone, cz) = cz.index();
3386 }
3387
3388 // Walk to split into regions
3389 const regionSplit cellRegion(mesh, isBlockedFace);
3390
3391 // Count number of cells per zone and per region
3392 labelList nCellsPerRegion(cellRegion.nRegions(), 0);
3393 labelList regionToZone(cellRegion.nRegions(), -2);
3394 labelList nCellsPerZone(czm.size()+1, 0);
3395 forAll(cellRegion, celli)
3396 {
3397 const label regioni = cellRegion[celli];
3398 const label zonei = cellToZone[celli];
3399
3400 // Zone for this region
3401 regionToZone[regioni] = zonei;
3402
3403 nCellsPerRegion[regioni]++;
3404 nCellsPerZone[zonei]++;
3405 }
3406 Pstream::listReduce(nCellsPerRegion, sumOp<label>());
3407 Pstream::listReduce(regionToZone, maxOp<label>());
3408 Pstream::listReduce(nCellsPerZone, sumOp<label>());
3409
3410
3411 // Mark small regions. Note that all processors have the same information
3412 // so should take the same decision. Alternative is to do master only
3413 // and scatter the decision.
3414 forAll(nCellsPerRegion, regioni)
3415 {
3416 const label zonei = regionToZone[regioni];
3417 label& nRegionCells = nCellsPerRegion[regioni];
3418
3419 if
3420 (
3421 nRegionCells < refineParams.minCellFraction()*nCellsPerZone[zonei]
3422 || nRegionCells < refineParams.nMinCells()
3423 )
3424 {
3425 Info<< "Deleting region " << regioni
3426 << " (size " << nRegionCells
3427 << ") of zone size " << nCellsPerZone[zonei]
3428 << endl;
3429
3430 // Mark region to be deleted. 0 size (= global) should never
3431 // occur.
3432 nRegionCells = 0;
3433 }
3434 }
3435
3436 DynamicList<label> cellsToRemove(mesh.nCells()/128);
3437 forAll(cellRegion, celli)
3438 {
3439 if (nCellsPerRegion[cellRegion[celli]] == 0)
3440 {
3441 cellsToRemove.append(celli);
3442 }
3443 }
3444 const label nTotCellsToRemove = returnReduce
3445 (
3446 cellsToRemove.size(),
3447 sumOp<label>()
3448 );
3449 if (nTotCellsToRemove > 0)
3450 {
3451 Info<< "Deleting " << nTotCellsToRemove
3452 << " cells in small regions" << endl;
3453
3454 removeCells cellRemover(mesh);
3455
3456 cellsToRemove.shrink();
3457 const labelList exposedFaces
3458 (
3459 cellRemover.getExposedFaces(cellsToRemove)
3460 );
3461 const labelList exposedPatch
3462 (
3463 UIndirectList<label>(ownPatch, exposedFaces)
3464 );
3465 (void)meshRefiner_.doRemoveCells
3466 (
3467 cellsToRemove,
3468 exposedFaces,
3469 exposedPatch,
3470 cellRemover
3471 );
3472 }
3473}
3474
3475
3477(
3478 const dictionary& refineDict,
3479 const refinementParameters& refineParams,
3480 const snapParameters& snapParams,
3481 const bool prepareForSnapping,
3482 const meshRefinement::FaceMergeType mergeType,
3483 const dictionary& motionDict
3484)
3485{
3486 addProfiling(refine, "snappyHexMesh::refine");
3487 Info<< nl
3488 << "Refinement phase" << nl
3489 << "----------------" << nl
3490 << endl;
3491
3492 const fvMesh& mesh = meshRefiner_.mesh();
3493
3494
3495 // Check that all the keep points are inside the mesh.
3496 refineParams.findCells(true, mesh, refineParams.locationsInMesh());
3497
3498 // Check that all the keep points are inside the mesh.
3499 if (dryRun_)
3500 {
3501 refineParams.findCells(true, mesh, refineParams.locationsOutsideMesh());
3502 }
3503
3504 // Estimate cell sizes
3505 if (dryRun_)
3506 {
3507 snappyVoxelMeshDriver voxelDriver
3508 (
3509 meshRefiner_,
3510 globalToMasterPatch_,
3511 globalToSlavePatch_
3512 );
3513 voxelDriver.doRefine(refineParams);
3514 }
3515
3516
3517 // Refine around feature edges
3518 featureEdgeRefine
3519 (
3520 refineParams,
3521 100, // maxIter
3522 0 // min cells to refine
3523 );
3524
3525
3526 // Initial automatic gap-level refinement: gaps smaller than the cell size
3527 if
3528 (
3529 max(meshRefiner_.surfaces().maxGapLevel()) > 0
3530 || max(meshRefiner_.shells().maxGapLevel()) > 0
3531 || max(meshRefiner_.surfaces().maxCurvatureLevel()) > 0
3532 )
3533 {
3534 // In case we use automatic gap level refinement do some pre-refinement
3535 // (fast) since it is so slow.
3536
3537 // Refine based on surface
3538 surfaceOnlyRefine
3539 (
3540 refineParams,
3541 20, // maxIter
3542 labelMax // no leak detection yet since all in same cell
3543 );
3544
3545 // Refine cells that contain a gap
3546 smallFeatureRefine
3547 (
3548 refineParams,
3549 100 // maxIter
3550 );
3551 }
3552
3553
3554 // Refine based on surface
3555 surfaceOnlyRefine
3556 (
3557 refineParams,
3558 100, // maxIter
3559 0 // Iteration to start leak detection and closure
3560 );
3561
3562 // Pass1 of automatic gap-level refinement: surface-intersected cells
3563 // in narrow gaps. Done early so we can remove the inside
3564 gapOnlyRefine
3565 (
3566 refineParams,
3567 100 // maxIter
3568 );
3569
3570 // Remove cells inbetween two surfaces
3571 surfaceProximityBlock
3572 (
3573 refineParams,
3574 1 //100 // maxIter
3575 );
3576
3577 // Remove cells (a certain distance) beyond surface intersections
3578 removeInsideCells
3579 (
3580 refineParams,
3581 1 // nBufferLayers
3582 );
3583
3584 // Pass2 of automatic gap-level refinement: all cells in gaps
3585 bigGapOnlyRefine
3586 (
3587 refineParams,
3588 true, // spreadGapSize
3589 100 // maxIter
3590 );
3591
3592 // Internal mesh refinement
3593 shellRefine
3594 (
3595 refineParams,
3596 100 // maxIter
3597 );
3598
3599 // Remove any extra cells from limitRegion with level -1, without
3600 // adding any buffer layer this time
3601 removeInsideCells
3602 (
3603 refineParams,
3604 0 // nBufferLayers
3605 );
3606
3607 // Refine any hexes with 5 or 6 faces refined to make smooth edges
3608 danglingCellRefine
3609 (
3610 refineParams,
3611 21, // 1 coarse face + 5 refined faces
3612 100 // maxIter
3613 );
3614 danglingCellRefine
3615 (
3616 refineParams,
3617 24, // 0 coarse faces + 6 refined faces
3618 100 // maxIter
3619 );
3620
3621 // Refine any cells with differing refinement level on either side
3622 refinementInterfaceRefine
3623 (
3624 refineParams,
3625 10 // maxIter
3626 );
3627
3628 // Directional shell refinement
3629 directionalShellRefine
3630 (
3631 refineParams,
3632 100 // maxIter
3633 );
3634
3635 // Block gaps (always, ignore surface leakLevel)
3636 if (refineParams.locationsOutsideMesh().size())
3637 {
3638 // For now: only check leaks on meshed surfaces. The problem is that
3639 // blockLeakFaces always generates baffles any not just faceZones ...
3640 const labelList unnamedSurfaces
3641 (
3643 (
3644 meshRefiner_.surfaces().surfZones()
3645 )
3646 );
3647 if (unnamedSurfaces.size() && !dryRun_)
3648 {
3649 meshRefiner_.blockLeakFaces
3650 (
3651 globalToMasterPatch_,
3652 globalToSlavePatch_,
3653 refineParams.locationsInMesh(),
3654 refineParams.zonesInMesh(),
3655 refineParams.locationsOutsideMesh(),
3656 unnamedSurfaces
3657 );
3658 }
3659 }
3660
3661 if
3662 (
3663 max(meshRefiner_.shells().nSmoothExpansion()) > 0
3664 || max(meshRefiner_.shells().nSmoothPosition()) > 0
3665 )
3666 {
3667 directionalSmooth(refineParams);
3668 }
3669
3670
3671 // Introduce baffles at surface intersections. Remove sections unreachable
3672 // from keepPoint.
3673 baffleAndSplitMesh
3674 (
3675 refineParams,
3676 snapParams,
3677 prepareForSnapping,
3678 motionDict
3679 );
3680
3681
3682 //- Commented out for now since causes zoning errors (sigsegv) on
3683 // case with faceZones. TBD.
3686 //boundaryRefinementInterfaceRefine
3687 //(
3688 // refineParams,
3689 // 10 // maxIter
3690 //);
3691
3692
3693 // Mesh is at its finest. Do optional zoning (cellZones and faceZones)
3694 wordPairHashTable zonesToFaceZone;
3695 zonify(refineParams, zonesToFaceZone);
3696
3697 // Create pairs of patches for faceZones
3698 {
3699 HashTable<Pair<word>> faceZoneToPatches(zonesToFaceZone.size());
3700
3701 // Note: zonesToFaceZone contains the same data on different
3702 // processors but in different order. We could sort the
3703 // contents but instead just loop in sortedToc order.
3704 List<Pair<word>> czs(zonesToFaceZone.sortedToc());
3705
3706 forAll(czs, i)
3707 {
3708 const Pair<word>& czNames = czs[i];
3709 const word& fzName = zonesToFaceZone[czNames];
3710
3711 const word& masterName = fzName;
3712 const word slaveName = czNames.second() + "_to_" + czNames.first();
3713 Pair<word> patches(masterName, slaveName);
3714 faceZoneToPatches.insert(fzName, patches);
3715 }
3716 addFaceZones(meshRefiner_, refineParams, faceZoneToPatches);
3717 }
3718
3719 // Pull baffles apart
3720 splitAndMergeBaffles
3721 (
3722 refineParams,
3723 snapParams,
3724 prepareForSnapping,
3725 motionDict
3726 );
3727
3728
3729 // Do something about cells with refined faces on the boundary
3730 if (prepareForSnapping)
3731 {
3732 const auto mt = meshRefiner_.meshType();
3733 if
3734 (
3735 meshRefiner_.mesh().faceZones().size()
3736 && (
3739 )
3740 )
3741 {
3742 erodeNonManifoldZoneFaces(refineParams);
3743 }
3744
3745 mergePatchFaces(mergeType, refineParams, motionDict);
3746 }
3747
3748
3749 if (refineParams.minCellFraction() > 0 || refineParams.nMinCells() > 0)
3750 {
3751 // Some small disconnected bits of mesh might remain since at
3752 // this point faceZones have not been converted into e.g. baffles.
3753 // We don't know whether e.g. the baffles are reset to be cyclicAMI
3754 // thus reconnecting. For now check if there are any particularly
3755 // small regions.
3756 deleteSmallRegions(refineParams);
3757 }
3758
3759
3760 if (!dryRun_ && Pstream::parRun())
3761 {
3762 Info<< nl
3763 << "Doing final balancing" << nl
3764 << "---------------------" << nl
3765 << endl;
3766
3767
3768 /*
3769 const bool hasBufferLayer
3770 (
3771 (meshRefiner_.meshType() == meshRefinement::CASTELLATEDBUFFERLAYER)
3772 || (meshRefiner_.meshType() == meshRefinement::CASTELLATEDBUFFERLAYER2)
3773 );
3774 */
3775 const bool hasBufferLayer = false;
3776
3777 labelList singleProcPoints;
3778
3779 /*
3780 if (hasBufferLayer)
3781 {
3782 //- Needed since buffer layer addition did not handle
3783 //- inter-processor extrusion. Fixed now. This code can be removed.
3784
3785 // Pick up all points on surfaces with specified buffer layers
3786 const labelList& surfIndex = meshRefiner_.surfaceIndex();
3787 const auto& addLayers = meshRefiner_.surfaces().addBufferLayers();
3788
3789 bitSet isSelected(mesh.nPoints());
3790 forAll(surfIndex, facei)
3791 {
3792 const label surfi = surfIndex[facei];
3793 if (surfIndex[facei] != -1)
3794 {
3795 const label globalRegioni =
3796 meshRefiner_.surfaces().globalRegion(surfi, 0);
3797 if (addLayers[globalRegioni])
3798 {
3799 isSelected.set(mesh.faces()[facei]);
3800 }
3801 }
3802 }
3803 syncTools::syncPointList
3804 (
3805 mesh,
3806 isSelected,
3807 orEqOp<unsigned int>(),
3808 0
3809 );
3810 singleProcPoints = isSelected.sortedToc();
3811 }
3812 */
3813
3814 // Do final balancing. Keep zoned faces on one processor since the
3815 // snap phase will convert them to baffles and this only works for
3816 // internal faces.
3817 meshRefiner_.balance
3818 (
3819 true, // keepZoneFaces
3820 hasBufferLayer, // keepBaffles
3821 singleProcPoints, // keepZonePoints
3822 scalarField(mesh.nCells(), 1), // cellWeights
3823 decomposer_,
3824 distributor_
3825 );
3826 }
3827}
3828
3829
3830// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
scalar range
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
A List with indirect addressing.
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
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
const labelListList & faceEdges() const
Return face-edge addressing.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
static void broadcasts(const int communicator, Type &value, Args &&... values)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A list that is sorted upon construction or when explicitly requested with the sort() method.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
T * first()
The first entry in the list.
Definition UILList.H:542
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const char *const componentNames[]
static const Form zero
static constexpr direction nComponents
Number of components in this vector space.
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp, const bitSet &orientation)
Per patch edge the pp faces (in global indices) using it.
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
labelList sortedToc() const
The indices of the on bits as a sorted labelList.
Definition bitSetI.H:441
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
A collection of cell labels.
Definition cellSet.H:50
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
Base class for writing coordSet(s) and tracks with fields.
double cpuTimeIncrement() const
Return CPU time [seconds] since last call to cpuTimeIncrement(), resetCpuTimeIncrement().
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
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
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
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label nTotalCells() const noexcept
Total global number of mesh cells.
Refinement of (split) hexes using polyTopoChange.
Definition hexRef8.H:63
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
label addFaceZone(const word &fzName, const word &masterPatch, const word &slavePatch, const surfaceZonesInfo::faceZoneType &fzType)
Add/lookup faceZone and update information. Return index of.
writeType
Enumeration for what to write. Used as a bit-pattern.
@ REMOVE
set value to -1 any face that was refined
debugType
Enumeration for what to debug. Used as a bit-pattern.
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
const fvMesh & mesh() const
Reference to mesh.
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
static writeType writeLevel()
Get/set write level.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label nNonProcessor() const
The number of patches before the first processor patch.
const labelList & patchID() const
Per boundary face label the patch index.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nInternalFaces() const noexcept
Number of internal faces.
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const vectorField & faceAreas() const
const cellList & cells() const
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
Simple container to keep together refinement specific information.
scalar planarAngle() const
Angle when two intersections are considered to be planar.
label maxLocalCells() const
Per processor max number of cells.
scalar curvature() const
Curvature.
scalar minCellFraction() const
When are disconnected regions small. Fraction of overall size.
label minRefineCells() const
When to stop refining.
label nMinCells() const
When are disconnected regions small. Absolute number of cells.
const pointField & locationsOutsideMesh() const
Optional points which are checked to be outside the mesh.
static labelList findCells(const bool checkInsideMesh, const polyMesh &, const pointField &locations)
Checks that cells are in mesh. Returns cells (or -1) they.
const pointField & locationsInMesh() const
Areas to keep.
dictionary getZoneInfo(const word &fzName, surfaceZonesInfo::faceZoneType &faceType) const
Get patchInfo and faceType for faceZone.
label maxCellUnbalance() const
Trigger cell count to start balancing.
scalar maxLoadUnbalance() const
Allowed load unbalance.
label maxGlobalCells() const
Total number of cells.
const wordList & zonesInMesh() const
Per area the zone name.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
const labelList & surfaces() const
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Given list of cells to remove, insert all the topology changes.
Definition removeCells.H:60
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Encapsulates queries for volume refinement ('refine all cells within shell').
const labelList & shells() const
Indices of surfaces that are shells.
Simple container to keep together snap specific information.
static void addFaceZones(meshRefinement &meshRefiner, const refinementParameters &refineParams, const HashTable< Pair< word > > &faceZoneToPatches)
Helper: add faceZones and patches.
void doRefine(const dictionary &refineDict, const refinementParameters &refineParams, const snapParameters &snapParams, const bool prepareForSnapping, const meshRefinement::FaceMergeType mergeType, const dictionary &motionDict)
Do all the refinement.
Equivalent of snappyRefineDriver but operating on a voxel mesh.
void doRefine(const refinementParameters &refineParams)
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName).
faceZoneType
What to do with faceZone faces.
static const Enum< faceZoneType > faceZoneTypeNames
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of named surfaces (surfaces with faceZoneName).
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
@ INSIDE
A location inside the volume.
Definition volumeType.H:65
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
static bool split(const std::string &line, std::string &key, std::string &val)
Definition cpuInfo.C:32
bool coupled
dynamicFvMesh & mesh
auto & name
const cellShapeList & cells
const labelIOList & zoneIDs
Definition correctPhi.H:59
const bitSet isBlockedFace(intersectedFaces())
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition ListListOps.C:62
List< T > createWithValue(const label len, const labelUList &locations, const T &val, const T &deflt=T())
Create a List filled with default values and various locations with another specified value.
@ CLAMP
Clamp value to the start/end value.
Definition tableBounds.H:65
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
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
HashTable< word, wordPair, Foam::Hash< wordPair > > wordPairHashTable
HashTable of wordPair.
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Vector< label > labelVector
Vector of labels.
Definition labelVector.H:47
constexpr label labelMin
Definition label.H:54
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
Smanip< std::ios_base::fmtflags > setf(std::ios_base::fmtflags flags)
Definition IOmanip.H:169
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
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
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
DynamicID< faceZoneMesh > faceZoneID
Foam::faceZoneID.
Definition ZoneIDs.H:43
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1, const label comm)
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
wordList patchNames(nPatches)
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
volScalarField & e
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Object access operator or list access operator (default is pass-through).
Definition UList.H:1120
Unit conversion functions.