Loading...
Searching...
No Matches
extrudeToRegionMesh.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 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
27Application
28 extrudeToRegionMesh
29
30Group
31 grpMeshGenerationUtilities
32
33Description
34 Extrude faceZones (internal or boundary faces) or faceSets (boundary faces
35 only) into a separate mesh (as a different region).
36
37 - used to e.g. extrude baffles (extrude internal faces) or create
38 liquid film regions.
39 - if extruding internal faces:
40 - create baffles in original mesh with mappedWall patches
41 - if extruding boundary faces:
42 - convert boundary faces to mappedWall patches
43 - extrude edges of faceZone as a <zone>_sidePatch
44 - extrude edges inbetween different faceZones as a
45 (nonuniformTransform)cyclic <zoneA>_<zoneB>
46 - extrudes into master direction (i.e. away from the owner cell
47 if flipMap is false)
48
49\verbatim
50
51Internal face extrusion
52-----------------------
53
54 +-------------+
55 | |
56 | |
57 +---AAAAAAA---+
58 | |
59 | |
60 +-------------+
61
62 AAA=faceZone to extrude.
63
64
65For the case of no flipMap the extrusion starts at owner and extrudes
66into the space of the neighbour:
67
68 +CCCCCCC+
69 | | <= extruded mesh
70 +BBBBBBB+
71
72 +-------------+
73 | |
74 | (neighbour) |
75 |___CCCCCCC___| <= original mesh (with 'baffles' added)
76 | BBBBBBB |
77 |(owner side) |
78 | |
79 +-------------+
80
81 BBB=mapped between owner on original mesh and new extrusion.
82 (zero offset)
83 CCC=mapped between neighbour on original mesh and new extrusion
84 (offset due to the thickness of the extruded mesh)
85
86For the case of flipMap the extrusion is the other way around: from the
87neighbour side into the owner side.
88
89
90Boundary face extrusion
91-----------------------
92
93 +--AAAAAAA--+
94 | |
95 | |
96 +-----------+
97
98 AAA=faceZone to extrude. E.g. slave side is owner side (no flipmap)
99
100becomes
101
102 +CCCCCCC+
103 | | <= extruded mesh
104 +BBBBBBB+
105
106 +--BBBBBBB--+
107 | | <= original mesh
108 | |
109 +-----------+
110
111 BBB=mapped between original mesh and new extrusion
112 CCC=polypatch
113
114
115Notes:
116 - when extruding cyclics with only one cell inbetween it does not
117 detect this as a cyclic since the face is the same face. It will
118 only work if the coupled edge extrudes a different face so if there
119 are more than 1 cell inbetween.
120
121\endverbatim
122
123\*---------------------------------------------------------------------------*/
124
125#include "argList.H"
126#include "fvMesh.H"
127#include "polyTopoChange.H"
128#include "OFstream.H"
129#include "meshTools.H"
130#include "mappedWallPolyPatch.H"
131#include "createShellMesh.H"
132#include "syncTools.H"
133#include "cyclicPolyPatch.H"
134#include "wedgePolyPatch.H"
136#include "extrudeModel.H"
137#include "globalIndex.H"
138#include "faceSet.H"
139
140#include "volFields.H"
141#include "surfaceFields.H"
142#include "pointFields.H"
143//#include "ReadFields.H"
144#include "fvMeshTools.H"
145#include "OBJstream.H"
146#include "PatchTools.H"
147#include "processorMeshes.H"
148
149using namespace Foam;
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153label findPatchID(const UList<polyPatch*>& newPatches, const word& name)
154{
155 forAll(newPatches, i)
156 {
157 if (newPatches[i]->name() == name)
158 {
159 return i;
160 }
161 }
162 return -1;
163}
164
165
166#ifdef FULLDEBUG
167void printPatches(Ostream& os, const UList<polyPatch*>& newPatches)
168{
169 for (const polyPatch* ppPtr : newPatches)
170 {
171 const polyPatch& pp = *ppPtr;
172 os << " name:" << pp.name() << " index:" << pp.index()
173 << " start:" << pp.start() << " size:" << pp.size()
174 << " type:" << pp.type() << nl;
175 }
176}
177#endif
178
179
180template<class PatchType>
181label addPatch
182(
184 const word& patchName,
185 DynamicList<polyPatch*>& newPatches
186)
187{
188 label patchi = findPatchID(newPatches, patchName);
189
190 if (patchi != -1)
191 {
192 if (isA<PatchType>(*newPatches[patchi]))
193 {
194 // Already there
195 return patchi;
196 }
197 else
198 {
200 << "Already have patch " << patchName
201 << " but of type " << newPatches[patchi]->type()
202 << exit(FatalError);
203 }
204 }
205
206
207 patchi = newPatches.size();
208
209 label startFacei = 0;
210 if (patchi > 0)
211 {
212 const polyPatch& pp = *newPatches.last();
213 startFacei = pp.start()+pp.size();
214 }
215
216 #ifdef FULLDEBUG
217 Pout<< "addPatch : starting newPatches:"
218 << " patch:" << patchi << " startFace:" << startFacei << nl;
219 printPatches(Pout, newPatches);
220 Pout<< "*** end of addPatch:" << endl;
221 #endif
222
223 newPatches.append
224 (
226 (
227 PatchType::typeName,
228 patchName,
229 0, // size
230 startFacei, // nFaces
231 patchi,
232 patches
233 ).ptr()
234 );
235
236 return patchi;
237}
238
239
240template<class PatchType>
241label addPatch
242(
244 const word& patchName,
245 const dictionary& dict,
246 DynamicList<polyPatch*>& newPatches
247)
248{
249 label patchi = findPatchID(newPatches, patchName);
250
251 if (patchi != -1)
252 {
253 if (isA<PatchType>(*newPatches[patchi]))
254 {
255 // Already there
256 return patchi;
257 }
258 else
259 {
261 << "Already have patch " << patchName
262 << " but of type " << newPatches[patchi]->type()
263 << exit(FatalError);
264 }
265 }
266
267
268 patchi = newPatches.size();
269
270 label startFacei = 0;
271 if (patchi > 0)
272 {
273 const polyPatch& pp = *newPatches.last();
274 startFacei = pp.start()+pp.size();
275 }
276
277
278 #ifdef FULLDEBUG
279 Pout<< "addPatch : starting newPatches:"
280 << " patch:" << patchi << " startFace:" << startFacei << nl;
281 printPatches(Pout, newPatches);
282 Pout<< "*** end of addPatch:" << endl;
283 #endif
284
285
286 dictionary patchDict(dict);
287 patchDict.set("type", PatchType::typeName);
288 patchDict.set("nFaces", 0);
289 patchDict.set("startFace", startFacei);
290
291 newPatches.append
292 (
294 (
295 patchName,
296 patchDict,
297 patchi,
298 patches
299 ).ptr()
300 );
301
302 return patchi;
303}
304
305
306// Remove zero-sized patches
307void deleteEmptyPatches(fvMesh& mesh)
308{
309 const polyBoundaryMesh& patches = mesh.boundaryMesh();
310
311 wordList masterNames;
312 if (Pstream::master())
313 {
314 masterNames = patches.names();
315 }
316 Pstream::broadcast(masterNames);
317
318
319 labelList oldToNew(patches.size(), -1);
320 label usedI = 0;
321 label notUsedI = patches.size();
322
323 // Add all the non-empty, non-processor patches
324 forAll(masterNames, masterI)
325 {
326 label patchi = patches.findPatchID(masterNames[masterI]);
327
328 if (patchi != -1)
329 {
330 if (isA<processorPolyPatch>(patches[patchi]))
331 {
332 // Similar named processor patch? Not 'possible'.
333 if (patches[patchi].empty())
334 {
335 Pout<< "Deleting processor patch " << patchi
336 << " name:" << patches[patchi].name()
337 << endl;
338 oldToNew[patchi] = --notUsedI;
339 }
340 else
341 {
342 oldToNew[patchi] = usedI++;
343 }
344 }
345 else
346 {
347 // Common patch.
348 if (returnReduceAnd(patches[patchi].empty()))
349 {
350 Pout<< "Deleting patch " << patchi
351 << " name:" << patches[patchi].name()
352 << endl;
353 oldToNew[patchi] = --notUsedI;
354 }
355 else
356 {
357 oldToNew[patchi] = usedI++;
358 }
359 }
360 }
361 }
362
363 // Add remaining patches at the end
364 forAll(patches, patchi)
365 {
366 if (oldToNew[patchi] == -1)
367 {
368 // Unique to this processor. Note: could check that these are
369 // only processor patches.
370 if (patches[patchi].empty())
371 {
372 Pout<< "Deleting processor patch " << patchi
373 << " name:" << patches[patchi].name()
374 << endl;
375 oldToNew[patchi] = --notUsedI;
376 }
377 else
378 {
379 oldToNew[patchi] = usedI++;
380 }
381 }
382 }
383
384 fvMeshTools::reorderPatches(mesh, oldToNew, usedI, true);
385}
386
387
388// Check zone either all internal or all external faces
389void checkZoneInside
390(
391 const polyMesh& mesh,
392 const wordList& zoneNames,
393 const labelList& zoneID,
394 const labelList& extrudeMeshFaces,
395 const boolList& isInternal
396)
397{
398 forAll(zoneNames, i)
399 {
400 if (isInternal[i])
401 {
402 Info<< "Zone " << zoneNames[i] << " has internal faces" << endl;
403 }
404 else
405 {
406 Info<< "Zone " << zoneNames[i] << " has boundary faces" << endl;
407 }
408 }
409
410 forAll(extrudeMeshFaces, i)
411 {
412 label facei = extrudeMeshFaces[i];
413 label zoneI = zoneID[i];
414 if (isInternal[zoneI] != mesh.isInternalFace(facei))
415 {
417 << "Zone " << zoneNames[zoneI]
418 << " is not consistently all internal or all boundary faces."
419 << " Face " << facei << " at " << mesh.faceCentres()[facei]
420 << " is the first occurrence."
421 << exit(FatalError);
422 }
423 }
424}
425
426
427// Calculate global pp faces per pp edge.
428labelListList globalEdgeFaces
429(
430 const polyMesh& mesh,
431 const globalIndex& globalFaces,
432 const primitiveFacePatch& pp,
433 const labelList& ppMeshEdges
434)
435{
436 // From mesh edge to global pp face labels.
437 labelListList globalEdgeFaces(ppMeshEdges.size());
438
439 const labelListList& edgeFaces = pp.edgeFaces();
440
441 forAll(edgeFaces, edgeI)
442 {
443 // Store pp face and processor as unique tag.
444 globalEdgeFaces[edgeI] = globalFaces.toGlobal(edgeFaces[edgeI]);
445 }
446
447 // Synchronise across coupled edges.
449 (
450 mesh,
451 ppMeshEdges,
452 globalEdgeFaces,
454 labelList() // null value
455 );
456
457 return globalEdgeFaces;
458}
459
460
461// Find a patch face that is not extruded. Return -1 if not found.
462label findUncoveredPatchFace
463(
464 const fvMesh& mesh,
465 const labelUIndList& extrudeMeshFaces, // mesh faces that are extruded
466 const label meshEdgeI // mesh edge
467)
468{
469 // Make set of extruded faces.
470 labelHashSet extrudeFaceSet(extrudeMeshFaces.size());
471 extrudeFaceSet.insert(extrudeMeshFaces);
472
473 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
474 const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
475 forAll(eFaces, i)
476 {
477 label facei = eFaces[i];
478 label patchi = pbm.whichPatch(facei);
479
480 if
481 (
482 patchi != -1
483 && !pbm[patchi].coupled()
484 && !extrudeFaceSet.found(facei)
485 )
486 {
487 return facei;
488 }
489 }
490
491 return -1;
492}
493
494
495// Same as findUncoveredPatchFace, except explicitly checks for cyclic faces
496label findUncoveredCyclicPatchFace
497(
498 const fvMesh& mesh,
499 const labelUIndList& extrudeMeshFaces, // mesh faces that are extruded
500 const label meshEdgeI // mesh edge
501)
502{
503 // Make set of extruded faces.
504 labelHashSet extrudeFaceSet(extrudeMeshFaces.size());
505 extrudeFaceSet.insert(extrudeMeshFaces);
506
507 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
508 const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
509 forAll(eFaces, i)
510 {
511 label facei = eFaces[i];
512 label patchi = pbm.whichPatch(facei);
513
514 if
515 (
516 patchi != -1
517 && isA<cyclicPolyPatch>(pbm[patchi])
518 && !extrudeFaceSet.found(facei)
519 )
520 {
521 return facei;
522 }
523 }
524
525 return -1;
526}
527
528
529// Calculate per edge min and max zone
530void calcEdgeMinMaxZone
531(
532 const fvMesh& mesh,
533 const primitiveFacePatch& extrudePatch,
534 const labelList& extrudeMeshEdges,
535 const labelList& zoneID,
536 const mapDistribute& extrudeEdgeFacesMap,
537 const labelListList& extrudeEdgeGlobalFaces,
538
539 labelList& minZoneID,
540 labelList& maxZoneID
541)
542{
543 // Get zoneIDs in extrudeEdgeGlobalFaces order
544 labelList mappedZoneID(zoneID);
545 extrudeEdgeFacesMap.distribute(mappedZoneID);
546
547 // Get min and max zone per edge
548 minZoneID.setSize(extrudeEdgeGlobalFaces.size(), labelMax);
549 maxZoneID.setSize(extrudeEdgeGlobalFaces.size(), labelMin);
550
551 forAll(extrudeEdgeGlobalFaces, edgeI)
552 {
553 const labelList& eFaces = extrudeEdgeGlobalFaces[edgeI];
554 if (eFaces.size())
555 {
556 forAll(eFaces, i)
557 {
558 label zoneI = mappedZoneID[eFaces[i]];
559 minZoneID[edgeI] = Foam::min(minZoneID[edgeI], zoneI);
560 maxZoneID[edgeI] = Foam::max(maxZoneID[edgeI], zoneI);
561 }
562 }
563 }
565 (
566 mesh,
567 extrudeMeshEdges,
568 minZoneID,
570 labelMax // null value
571 );
573 (
574 mesh,
575 extrudeMeshEdges,
576 maxZoneID,
578 labelMin // null value
579 );
580}
581
582
583// Count the number of faces in patches that need to be created. Calculates:
584// zoneSidePatch[zoneI] : the number of side faces to be created
585// zoneZonePatch[zoneA,zoneB] : the number of faces inbetween zoneA and B
586// Since this only counts we're not taking the processor patches into
587// account.
588void countExtrudePatches
589(
590 const fvMesh& mesh,
591 const label nZones,
592 const primitiveFacePatch& extrudePatch,
593 const labelList& extrudeMeshFaces,
594 const labelList& extrudeMeshEdges,
595
596 const labelListList& extrudeEdgeGlobalFaces,
597 const labelList& minZoneID,
598 const labelList& maxZoneID,
599
600 labelList& zoneSidePatch,
601 labelList& zoneZonePatch
602)
603{
604 // Check on master edge for use of zones. Since we only want to know
605 // whether they are being used at all no need to accurately count on slave
606 // edge as well. Just add all together at the end of this routine so it
607 // gets detected at least.
608
609 forAll(extrudePatch.edgeFaces(), edgeI)
610 {
611 const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
612
613 if (eFaces.size() == 2)
614 {
615 // Internal edge - check if inbetween different zones.
616 if (minZoneID[edgeI] != maxZoneID[edgeI])
617 {
618 zoneZonePatch[minZoneID[edgeI]*nZones+maxZoneID[edgeI]]++;
619 }
620 }
621 else if
622 (
623 eFaces.size() == 1
624 && extrudeEdgeGlobalFaces[edgeI].size() == 2
625 )
626 {
627 // Coupled edge - check if inbetween different zones.
628 if (minZoneID[edgeI] != maxZoneID[edgeI])
629 {
630 const edge& e = extrudePatch.edges()[edgeI];
631 const pointField& pts = extrudePatch.localPoints();
633 << "Edge " << edgeI
634 << "at " << pts[e[0]] << pts[e[1]]
635 << " is a coupled edge and inbetween two different zones "
636 << minZoneID[edgeI] << " and " << maxZoneID[edgeI] << endl
637 << " This is currently not supported." << endl;
638
639 zoneZonePatch[minZoneID[edgeI]*nZones+maxZoneID[edgeI]]++;
640 }
641 }
642 else
643 {
644 // One or more than two edge-faces.
645 // Check whether we are on a mesh edge with external patches. If
646 // so choose any uncovered one. If none found put face in
647 // undetermined zone 'side' patch
648
649 label facei = findUncoveredPatchFace
650 (
651 mesh,
652 labelUIndList(extrudeMeshFaces, eFaces),
653 extrudeMeshEdges[edgeI]
654 );
655
656 if (facei == -1)
657 {
658 zoneSidePatch[minZoneID[edgeI]]++;
659 }
660 }
661 }
662 // Synchronise decision. Actual numbers are not important, just make
663 // sure that they're > 0 on all processors.
664 Pstream::listReduce(zoneSidePatch, sumOp<label>());
665 Pstream::listReduce(zoneZonePatch, sumOp<label>());
666}
667
668
669void addCouplingPatches
670(
671 const fvMesh& mesh,
672 const word& regionName,
673 const word& shellRegionName,
674 const wordList& zoneNames,
675 const wordList& zoneShadowNames,
676 const boolList& isInternal,
677 const labelList& zoneIDs,
678
679 DynamicList<polyPatch*>& newPatches,
680 labelList& interRegionTopPatch,
681 labelList& interRegionBottomPatch
682)
683{
684 Pout<< "Adding coupling patches:" << nl << nl
685 << "patchID\tpatch\ttype" << nl
686 << "-------\t-----\t----"
687 << endl;
688
689 interRegionTopPatch.setSize(zoneNames.size(), -1);
690 interRegionBottomPatch.setSize(zoneNames.size(), -1);
691
692 label nOldPatches = newPatches.size();
693 forAll(zoneNames, zoneI)
694 {
695 word interName
696 (
698 +"_to_"
699 +shellRegionName
700 +'_'
701 +zoneNames[zoneI]
702 );
703
704 if (isInternal[zoneI])
705 {
706 interRegionTopPatch[zoneI] = addPatch<mappedWallPolyPatch>
707 (
708 mesh.boundaryMesh(),
709 interName + "_top",
710 newPatches
711 );
712 Pout<< interRegionTopPatch[zoneI]
713 << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
714 << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
715 << nl;
716
717 interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
718 (
719 mesh.boundaryMesh(),
720 interName + "_bottom",
721 newPatches
722 );
723 Pout<< interRegionBottomPatch[zoneI]
724 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
725 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
726 << nl;
727 }
728 else if (zoneShadowNames.size() == 0)
729 {
730 interRegionTopPatch[zoneI] = addPatch<polyPatch>
731 (
732 mesh.boundaryMesh(),
733 zoneNames[zoneI] + "_top",
734 newPatches
735 );
736 Pout<< interRegionTopPatch[zoneI]
737 << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
738 << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
739 << nl;
740
741 interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
742 (
743 mesh.boundaryMesh(),
744 interName,
745 newPatches
746 );
747 Pout<< interRegionBottomPatch[zoneI]
748 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
749 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
750 << nl;
751 }
752 else //patch using shadow face zones.
753 {
754 interRegionTopPatch[zoneI] = addPatch<mappedWallPolyPatch>
755 (
756 mesh.boundaryMesh(),
757 zoneShadowNames[zoneI] + "_top",
758 newPatches
759 );
760 Pout<< interRegionTopPatch[zoneI]
761 << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
762 << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
763 << nl;
764
765 interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
766 (
767 mesh.boundaryMesh(),
768 interName,
769 newPatches
770 );
771 Pout<< interRegionBottomPatch[zoneI]
772 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
773 << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
774 << nl;
775 }
776 }
777 Pout<< "Added " << newPatches.size()-nOldPatches
778 << " inter-region patches." << nl
779 << endl;
780}
781
782
783// Sets sidePatch[edgeI] to interprocessor or cyclic patch. Adds any
784// coupled patches if necessary.
785void addCoupledPatches
786(
787 const fvMesh& mesh,
788 const primitiveFacePatch& extrudePatch,
789 const labelList& extrudeMeshFaces,
790 const labelList& extrudeMeshEdges,
791 const mapDistribute& extrudeEdgeFacesMap,
792 const labelListList& extrudeEdgeGlobalFaces,
793
794 labelList& sidePatchID,
795 DynamicList<polyPatch*>& newPatches
796)
797{
798 // Calculate opposite processor for coupled edges (only if shared by
799 // two procs). Note: could have saved original globalEdgeFaces structure.
800
801 // Get procID in extrudeEdgeGlobalFaces order
802 labelList procID(extrudeEdgeGlobalFaces.size(), Pstream::myProcNo());
803 extrudeEdgeFacesMap.distribute(procID);
804
805 labelList minProcID(extrudeEdgeGlobalFaces.size(), labelMax);
806 labelList maxProcID(extrudeEdgeGlobalFaces.size(), labelMin);
807
808 forAll(extrudeEdgeGlobalFaces, edgeI)
809 {
810 const labelList& eFaces = extrudeEdgeGlobalFaces[edgeI];
811 if (eFaces.size())
812 {
813 forAll(eFaces, i)
814 {
815 label proci = procID[eFaces[i]];
816 minProcID[edgeI] = Foam::min(minProcID[edgeI], proci);
817 maxProcID[edgeI] = Foam::max(maxProcID[edgeI], proci);
818 }
819 }
820 }
822 (
823 mesh,
824 extrudeMeshEdges,
825 minProcID,
827 labelMax // null value
828 );
830 (
831 mesh,
832 extrudeMeshEdges,
833 maxProcID,
835 labelMin // null value
836 );
837
838 Pout<< "Adding processor or cyclic patches:" << nl << nl
839 << "patchID\tpatch" << nl
840 << "-------\t-----"
841 << endl;
842
843 label nOldPatches = newPatches.size();
844
845 sidePatchID.setSize(extrudePatch.edgeFaces().size(), -1);
846 forAll(extrudePatch.edgeFaces(), edgeI)
847 {
848 const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
849
850 if
851 (
852 eFaces.size() == 1
853 && extrudeEdgeGlobalFaces[edgeI].size() == 2
854 )
855 {
856 // coupled boundary edge. Find matching patch.
857 label nbrProci = minProcID[edgeI];
858 if (nbrProci == Pstream::myProcNo())
859 {
860 nbrProci = maxProcID[edgeI];
861 }
862
863
864 if (nbrProci == Pstream::myProcNo())
865 {
866 // Cyclic patch since both procs the same. This cyclic should
867 // already exist in newPatches so no adding necessary.
868
869 label facei = findUncoveredCyclicPatchFace
870 (
871 mesh,
872 labelUIndList(extrudeMeshFaces, eFaces),
873 extrudeMeshEdges[edgeI]
874 );
875
876 if (facei != -1)
877 {
878 const polyBoundaryMesh& patches = mesh.boundaryMesh();
879
880 label newPatchi = findPatchID
881 (
882 newPatches,
883 patches[patches.whichPatch(facei)].name()
884 );
885
886 sidePatchID[edgeI] = newPatchi;
887 }
888 else
889 {
891 << "Unable to determine coupled patch addressing"
892 << abort(FatalError);
893 }
894 }
895 else
896 {
897 // Processor patch
898 word name
899 (
901 );
902
903 sidePatchID[edgeI] = findPatchID(newPatches, name);
904
905 if (sidePatchID[edgeI] == -1)
906 {
907 dictionary patchDict;
908 patchDict.add("myProcNo", Pstream::myProcNo());
909 patchDict.add("neighbProcNo", nbrProci);
910
911 sidePatchID[edgeI] = addPatch<processorPolyPatch>
912 (
913 mesh.boundaryMesh(),
914 name,
915 patchDict,
916 newPatches
917 );
918
919 Pout<< sidePatchID[edgeI] << '\t' << name
920 << nl;
921 }
922 }
923 }
924 }
925 Pout<< "Added " << newPatches.size()-nOldPatches
926 << " coupled patches." << nl
927 << endl;
928}
929
930
931void addZoneSidePatches
932(
933 const fvMesh& mesh,
934 const wordList& zoneNames,
935 const word& oneDPolyPatchType,
936
937 DynamicList<polyPatch*>& newPatches,
938 labelList& zoneSidePatch
939)
940{
941 Pout<< "Adding patches for sides on zones:" << nl << nl
942 << "patchID\tpatch" << nl
943 << "-------\t-----"
944 << endl;
945
946 label nOldPatches = newPatches.size();
947
948 forAll(zoneSidePatch, zoneI)
949 {
950 if (!oneDPolyPatchType.empty())
951 {
952 // Reuse single empty patch.
953 word patchName;
954 if (oneDPolyPatchType == "empty")
955 {
956 patchName = "oneDEmptyPatch";
957 zoneSidePatch[zoneI] = addPatch<emptyPolyPatch>
958 (
959 mesh.boundaryMesh(),
960 patchName,
961 newPatches
962 );
963 }
964 else if (oneDPolyPatchType == "wedge")
965 {
966 patchName = "oneDWedgePatch";
967 zoneSidePatch[zoneI] = addPatch<wedgePolyPatch>
968 (
969 mesh.boundaryMesh(),
970 patchName,
971 newPatches
972 );
973 }
974 else
975 {
977 << "Type " << oneDPolyPatchType << " does not exist "
978 << exit(FatalError);
979 }
980
981 Pout<< zoneSidePatch[zoneI] << '\t' << patchName << nl;
982 }
983 else if (zoneSidePatch[zoneI] > 0)
984 {
985 word patchName = zoneNames[zoneI] + "_" + "side";
986
987 zoneSidePatch[zoneI] = addPatch<polyPatch>
988 (
989 mesh.boundaryMesh(),
990 patchName,
991 newPatches
992 );
993
994 Pout<< zoneSidePatch[zoneI] << '\t' << patchName << nl;
995 }
996 }
997 Pout<< "Added " << newPatches.size()-nOldPatches << " zone-side patches."
998 << nl << endl;
999}
1000
1001
1002void addInterZonePatches
1003(
1004 const fvMesh& mesh,
1005 const wordList& zoneNames,
1006 const bool oneD,
1007
1008 labelList& zoneZonePatch_min,
1009 labelList& zoneZonePatch_max,
1010 DynamicList<polyPatch*>& newPatches
1011)
1012{
1013 Pout<< "Adding inter-zone patches:" << nl << nl
1014 << "patchID\tpatch" << nl
1015 << "-------\t-----"
1016 << endl;
1017
1018 dictionary transformDict;
1019 transformDict.add
1020 (
1021 "transform",
1023 );
1024
1025 label nOldPatches = newPatches.size();
1026
1027 if (!oneD)
1028 {
1029 forAll(zoneZonePatch_min, minZone)
1030 {
1031 for (label maxZone = minZone; maxZone < zoneNames.size(); maxZone++)
1032 {
1033 label index = minZone*zoneNames.size()+maxZone;
1034
1035 if (zoneZonePatch_min[index] > 0)
1036 {
1037 word minToMax =
1038 zoneNames[minZone]
1039 + "_to_"
1040 + zoneNames[maxZone];
1041 word maxToMin =
1042 zoneNames[maxZone]
1043 + "_to_"
1044 + zoneNames[minZone];
1045
1046 {
1047 transformDict.set("neighbourPatch", maxToMin);
1048 zoneZonePatch_min[index] =
1049 addPatch<nonuniformTransformCyclicPolyPatch>
1050 (
1051 mesh.boundaryMesh(),
1052 minToMax,
1053 transformDict,
1054 newPatches
1055 );
1056 Pout<< zoneZonePatch_min[index] << '\t' << minToMax
1057 << nl;
1058 }
1059 {
1060 transformDict.set("neighbourPatch", minToMax);
1061 zoneZonePatch_max[index] =
1062 addPatch<nonuniformTransformCyclicPolyPatch>
1063 (
1064 mesh.boundaryMesh(),
1065 maxToMin,
1066 transformDict,
1067 newPatches
1068 );
1069 Pout<< zoneZonePatch_max[index] << '\t' << maxToMin
1070 << nl;
1071 }
1072
1073 }
1074 }
1075 }
1076 }
1077 Pout<< "Added " << newPatches.size()-nOldPatches << " inter-zone patches."
1078 << nl << endl;
1079}
1080
1081
1082tmp<pointField> calcOffset
1083(
1084 const primitiveFacePatch& extrudePatch,
1085 const createShellMesh& extruder,
1086 const polyPatch& pp
1087)
1088{
1089 vectorField::subField fc = pp.faceCentres();
1090
1091 auto toffsets = tmp<pointField>::New(fc.size());
1092 auto& offsets = toffsets.ref();
1093
1094 forAll(fc, i)
1095 {
1096 label meshFacei = pp.start()+i;
1097 label patchFacei = mag(extruder.faceToFaceMap()[meshFacei])-1;
1098 point patchFc = extrudePatch[patchFacei].centre
1099 (
1100 extrudePatch.points()
1101 );
1102 offsets[i] = patchFc - fc[i];
1103 }
1104 return toffsets;
1105}
1106
1107
1108void setCouplingInfo
1109(
1110 fvMesh& mesh,
1111 const labelList& zoneToPatch,
1112 const word& sampleRegion,
1114 const List<pointField>& offsets,
1115 const List<boundBox>& bbs
1116)
1117{
1118 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1119
1120 List<polyPatch*> newPatches
1121 (
1122 patches.size(),
1123 static_cast<polyPatch*>(nullptr)
1124 );
1125
1126 forAll(zoneToPatch, zoneI)
1127 {
1128 const label patchi = zoneToPatch[zoneI];
1129
1130 if (patchi != -1)
1131 {
1132 const polyPatch& pp = patches[patchi];
1133
1135 {
1136 const vector avgOffset = gAverage(offsets[zoneI]);
1137 const scalar mergeSqrDist =
1138 gMax(magSqr(offsets[zoneI]-avgOffset));
1139
1140 // Create with uniform offset initially
1141 auto mappedPtr = autoPtr<mappedWallPolyPatch>::New
1142 (
1143 pp.name(),
1144 pp.size(),
1145 pp.start(),
1146 patchi,
1147 sampleRegion, // sampleRegion
1148 mode, // sampleMode
1149 pp.name(), // samplePatch
1150
1151 avgOffset, // uniform offset
1152 patches
1153 );
1154
1155 Info<< "Adding on " << mesh.name() << " coupling patch "
1156 << pp.name() << " with ";
1157
1158 // Verify uniformity of offset
1159 // (same check as blockMesh geom merge)
1160 if (mergeSqrDist < magSqr(10*SMALL*bbs[zoneI].span()))
1161 {
1162 Info<< "uniform offset " << avgOffset << endl;
1163 }
1164 else
1165 {
1166 Info<< "non-uniform offset" << endl;
1167
1168 (*mappedPtr).setOffset(offsets[zoneI]);
1169 }
1170
1171 newPatches[patchi] = mappedPtr.release();
1172 }
1173 }
1174 }
1175
1176 forAll(newPatches, patchi)
1177 {
1178 if (!newPatches[patchi])
1179 {
1180 newPatches[patchi] = patches[patchi].clone(patches).ptr();
1181 //newPatches[patchi].index() = patchi;
1182 }
1183 }
1184
1185
1186 #ifdef FULLDEBUG
1187 Pout<< "*** setCouplingInfo addFvPAtches:" << nl;
1188 printPatches(Pout, newPatches);
1189 Pout<< "*** setCouplingInfo end of addFvPAtches:" << endl;
1190 #endif
1191
1192 mesh.removeFvBoundary();
1193 mesh.addFvPatches(newPatches, true);
1194}
1195
1196
1197// Extrude and write geometric properties
1198void extrudeGeometricProperties
1199(
1200 const polyMesh& mesh,
1201 const primitiveFacePatch& extrudePatch,
1202 const createShellMesh& extruder,
1203 const polyMesh& regionMesh,
1204 const extrudeModel& model
1205)
1206{
1207 const pointIOField patchFaceCentres
1208 (
1209 IOobject
1210 (
1211 "patchFaceCentres",
1212 mesh.pointsInstance(),
1214 mesh,
1216 )
1217 );
1218
1219 const pointIOField patchEdgeCentres
1220 (
1221 IOobject
1222 (
1223 "patchEdgeCentres",
1224 mesh.pointsInstance(),
1226 mesh,
1228 )
1229 );
1230
1231 //forAll(extrudePatch.edges(), edgeI)
1232 //{
1233 // const edge& e = extrudePatch.edges()[edgeI];
1234 // Pout<< "Edge:" << e.centre(extrudePatch.localPoints()) << nl
1235 // << "read:" << patchEdgeCentres[edgeI]
1236 // << endl;
1237 //}
1238
1239
1240 // Determine edge normals on original patch
1241 labelList patchEdges;
1242 labelList coupledEdges;
1243 bitSet sameEdgeOrientation;
1245 (
1246 extrudePatch,
1247 mesh.globalData().coupledPatch(),
1248 patchEdges,
1249 coupledEdges,
1250 sameEdgeOrientation
1251 );
1252
1253 pointField patchEdgeNormals
1254 (
1256 (
1257 mesh,
1258 extrudePatch,
1259 patchEdges,
1260 coupledEdges
1261 )
1262 );
1263
1264
1265 pointIOField faceCentres
1266 (
1267 IOobject
1268 (
1269 "faceCentres",
1270 regionMesh.pointsInstance(),
1272 regionMesh,
1276 ),
1277 regionMesh.nFaces()
1278 );
1279
1280
1281 // Work out layers. Guaranteed in columns so no fancy parallel bits.
1282
1283
1284 forAll(extruder.faceToFaceMap(), facei)
1285 {
1286 if (extruder.faceToFaceMap()[facei] != 0)
1287 {
1288 // 'horizontal' face
1289 label patchFacei = mag(extruder.faceToFaceMap()[facei])-1;
1290
1291 label celli = regionMesh.faceOwner()[facei];
1292 if (regionMesh.isInternalFace(facei))
1293 {
1294 celli = Foam::max(celli, regionMesh.faceNeighbour()[facei]);
1295 }
1296
1297 // Calculate layer from cell numbering (see createShellMesh)
1298 label layerI = (celli % model.nLayers());
1299
1300 if
1301 (
1302 !regionMesh.isInternalFace(facei)
1303 && extruder.faceToFaceMap()[facei] > 0
1304 )
1305 {
1306 // Top face
1307 layerI++;
1308 }
1309
1310
1311 // Recalculate based on extrusion model
1312 faceCentres[facei] = model
1313 (
1314 patchFaceCentres[patchFacei],
1315 extrudePatch.faceNormals()[patchFacei],
1316 layerI
1317 );
1318 }
1319 else
1320 {
1321 // 'vertical face
1322 label patchEdgeI = extruder.faceToEdgeMap()[facei];
1323 label layerI =
1324 (
1325 regionMesh.faceOwner()[facei]
1326 % model.nLayers()
1327 );
1328
1329 // Extrude patch edge centre to this layer
1330 point pt0 = model
1331 (
1332 patchEdgeCentres[patchEdgeI],
1333 patchEdgeNormals[patchEdgeI],
1334 layerI
1335 );
1336 // Extrude patch edge centre to next layer
1337 point pt1 = model
1338 (
1339 patchEdgeCentres[patchEdgeI],
1340 patchEdgeNormals[patchEdgeI],
1341 layerI+1
1342 );
1343
1344 // Interpolate
1345 faceCentres[facei] = 0.5*(pt0+pt1);
1346 }
1347 }
1348
1349 pointIOField cellCentres
1350 (
1351 IOobject
1352 (
1353 "cellCentres",
1354 regionMesh.pointsInstance(),
1356 regionMesh,
1360 ),
1361 regionMesh.nCells()
1362 );
1363
1364 forAll(extruder.cellToFaceMap(), celli)
1365 {
1366 label patchFacei = extruder.cellToFaceMap()[celli];
1367
1368 // Calculate layer from cell numbering (see createShellMesh)
1369 label layerI = (celli % model.nLayers());
1370
1371 // Recalculate based on extrusion model
1372 point pt0 = model
1373 (
1374 patchFaceCentres[patchFacei],
1375 extrudePatch.faceNormals()[patchFacei],
1376 layerI
1377 );
1378 point pt1 = model
1379 (
1380 patchFaceCentres[patchFacei],
1381 extrudePatch.faceNormals()[patchFacei],
1382 layerI+1
1383 );
1384
1385 // Interpolate
1386 cellCentres[celli] = 0.5*(pt0+pt1);
1387 }
1388
1389
1390 // Bit of checking
1391 if (false)
1392 {
1393 OBJstream faceStr(regionMesh.time().path()/"faceCentres.obj");
1394 OBJstream cellStr(regionMesh.time().path()/"cellCentres.obj");
1395
1396 forAll(faceCentres, facei)
1397 {
1398 Pout<< "Model :" << faceCentres[facei] << endl
1399 << "regionMesh:" << regionMesh.faceCentres()[facei] << endl;
1400 faceStr.writeLine
1401 (
1402 faceCentres[facei],
1403 regionMesh.faceCentres()[facei]
1404 );
1405 }
1406 forAll(cellCentres, celli)
1407 {
1408 Pout<< "Model :" << cellCentres[celli] << endl
1409 << "regionMesh:" << regionMesh.cellCentres()[celli] << endl;
1410 cellStr.writeLine
1411 (
1412 cellCentres[celli],
1413 regionMesh.cellCentres()[celli]
1414 );
1415 }
1416 }
1417
1418
1419
1420 Info<< "Writing geometric properties for mesh " << regionMesh.name()
1421 << " to " << regionMesh.pointsInstance() << nl
1422 << endl;
1423
1424 bool ok = faceCentres.write() && cellCentres.write();
1425
1426 if (!ok)
1427 {
1429 << "Failed writing " << faceCentres.objectPath()
1430 << " and " << cellCentres.objectPath()
1431 << exit(FatalError);
1432 }
1433}
1434
1435
1436int main(int argc, char *argv[])
1437{
1439 (
1440 "Create region mesh by extruding a faceZone or faceSet"
1441 );
1442
1443 #include "addRegionOption.H"
1444 #include "addOverwriteOption.H"
1445
1447 (
1448 "dict", "file", "Alternative extrudeToRegionMeshDict"
1449 );
1450
1451 #include "setRootCase.H"
1452 #include "createTime.H"
1453 #include "createNamedMesh.H"
1454
1455 if (mesh.boundaryMesh().checkParallelSync(true))
1456 {
1457 List<wordList> allNames(Pstream::nProcs());
1458 allNames[Pstream::myProcNo()] = mesh.boundaryMesh().names();
1459 Pstream::allGatherList(allNames);
1460
1462 << "Patches are not synchronised on all processors."
1463 << " Per processor patches " << allNames
1464 << exit(FatalError);
1465 }
1466
1467
1468 const word oldInstance = mesh.pointsInstance();
1469 const bool overwrite = args.found("overwrite");
1470
1471 const word dictName("extrudeToRegionMeshDict");
1472
1474
1476
1477 // Point generator
1479
1480 // Region
1481 const word shellRegionName(dict.get<word>("region"));
1482
1483 // Faces to extrude - either faceZones or faceSets (boundary faces only)
1484 wordList zoneNames;
1485 wordList zoneShadowNames;
1486
1487 const bool hasZones = dict.found("faceZones");
1488 if (hasZones)
1489 {
1490 dict.readEntry("faceZones", zoneNames);
1491 dict.readIfPresent("faceZonesShadow", zoneShadowNames);
1492
1493 // Check
1494 if (dict.found("faceSets"))
1495 {
1497 << "Please supply faces to extrude either through 'faceZones'"
1498 << " or 'faceSets' entry. Found both."
1499 << exit(FatalIOError);
1500 }
1501 }
1502 else
1503 {
1504 dict.readEntry("faceSets", zoneNames);
1505 dict.readIfPresent("faceSetsShadow", zoneShadowNames);
1506 }
1507
1508
1509 mappedPatchBase::sampleMode sampleMode =
1510 mappedPatchBase::sampleModeNames_.get("sampleMode", dict);
1511
1512 const bool oneD(dict.get<bool>("oneD"));
1513 bool oneDNonManifoldEdges(false);
1514 word oneDPatchType(emptyPolyPatch::typeName);
1515 if (oneD)
1516 {
1517 oneDNonManifoldEdges = dict.getOrDefault("nonManifold", false);
1518 oneDPatchType = dict.get<word>("oneDPolyPatchType");
1519 }
1520
1521 const bool adaptMesh(dict.get<bool>("adaptMesh"));
1522
1523 const bool addSidePatches(dict.getOrDefault<bool>("addSidePatches", true));
1524
1525
1526
1527 if (hasZones)
1528 {
1529 Info<< "Extruding zones " << zoneNames
1530 << " on mesh " << regionName
1531 << " into shell mesh " << shellRegionName
1532 << endl;
1533 }
1534 else
1535 {
1536 Info<< "Extruding faceSets " << zoneNames
1537 << " on mesh " << regionName
1538 << " into shell mesh " << shellRegionName
1539 << endl;
1540 }
1541
1542 if (shellRegionName == regionName)
1543 {
1545 << "Cannot extrude into same region as mesh." << endl
1546 << "Mesh region : " << regionName << endl
1547 << "Shell region : " << shellRegionName
1548 << exit(FatalIOError);
1549 }
1550
1551
1552 if (oneD)
1553 {
1554 if (oneDNonManifoldEdges)
1555 {
1556 Info<< "Extruding as 1D columns with sides in patch type "
1557 << oneDPatchType
1558 << " and connected points (except on non-manifold areas)."
1559 << endl;
1560 }
1561 else
1562 {
1563 Info<< "Extruding as 1D columns with sides in patch type "
1564 << oneDPatchType
1565 << " and duplicated points (overlapping volumes)."
1566 << endl;
1567 }
1568 }
1569
1570
1571 if (addSidePatches && zoneNames.size() > 1)
1572 {
1573 Info<< "Extruding edges on more than one faceZone into boundary faces"
1574 << endl;
1575 }
1576 else
1577 {
1578 Info<< "Extruding internal edges into internal faces" << endl;
1579 }
1580
1581
1583 //IOobjectList objects(mesh, runTime.timeName());
1584 //
1586 //
1587 //PtrList<volScalarField> vsFlds;
1588 //ReadFields(mesh, objects, vsFlds);
1589 //
1590 //PtrList<volVectorField> vvFlds;
1591 //ReadFields(mesh, objects, vvFlds);
1592 //
1593 //PtrList<volSphericalTensorField> vstFlds;
1594 //ReadFields(mesh, objects, vstFlds);
1595 //
1596 //PtrList<volSymmTensorField> vsymtFlds;
1597 //ReadFields(mesh, objects, vsymtFlds);
1598 //
1599 //PtrList<volTensorField> vtFlds;
1600 //ReadFields(mesh, objects, vtFlds);
1601 //
1603 //
1604 //PtrList<surfaceScalarField> ssFlds;
1605 //ReadFields(mesh, objects, ssFlds);
1606 //
1607 //PtrList<surfaceVectorField> svFlds;
1608 //ReadFields(mesh, objects, svFlds);
1609 //
1610 //PtrList<surfaceSphericalTensorField> sstFlds;
1611 //ReadFields(mesh, objects, sstFlds);
1612 //
1613 //PtrList<surfaceSymmTensorField> ssymtFlds;
1614 //ReadFields(mesh, objects, ssymtFlds);
1615 //
1616 //PtrList<surfaceTensorField> stFlds;
1617 //ReadFields(mesh, objects, stFlds);
1618 //
1620 //
1621 //PtrList<pointScalarField> psFlds;
1622 //ReadFields(pointMesh::New(mesh), objects, psFlds);
1623 //
1624 //PtrList<pointVectorField> pvFlds;
1625 //ReadFields(pointMesh::New(mesh), objects, pvFlds);
1626
1627
1628
1629 // Create dummy fv* files
1630 fvMeshTools::createDummyFvMeshFiles(mesh, shellRegionName, true);
1631
1632
1633 word meshInstance;
1634 if (!overwrite)
1635 {
1636 ++runTime;
1637 meshInstance = runTime.timeName();
1638 }
1639 else
1640 {
1641 meshInstance = oldInstance;
1642 }
1643 Info<< "Writing meshes to " << meshInstance << nl << endl;
1644
1645
1646 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1647
1648
1649 // Extract faces to extrude
1650 // ~~~~~~~~~~~~~~~~~~~~~~~~
1651 // Note: zoneID are regions of extrusion. They are not mesh.faceZones
1652 // indices.
1653
1654 // From extrude zone to mesh zone (or -1 if extruding faceSets)
1655 labelList meshZoneID;
1656 // Per extrude zone whether contains internal or external faces
1657 boolList isInternal(zoneNames.size(), false);
1658
1659 labelList extrudeMeshFaces;
1660 faceList zoneFaces;
1661 labelList zoneID;
1662 boolList zoneFlipMap;
1663 // Shadow
1664 labelList zoneShadowIDs; // from extrude shadow zone to mesh zone
1665 labelList extrudeMeshShadowFaces;
1666 boolList zoneShadowFlipMap;
1667 labelList zoneShadowID;
1668
1669 if (hasZones)
1670 {
1671 const faceZoneMesh& faceZones = mesh.faceZones();
1672
1673 meshZoneID.setSize(zoneNames.size());
1674 forAll(zoneNames, i)
1675 {
1676 meshZoneID[i] = faceZones.findZoneID(zoneNames[i]);
1677 if (meshZoneID[i] == -1)
1678 {
1680 << "Cannot find zone " << zoneNames[i] << endl
1681 << "Valid zones are " << faceZones.names()
1682 << exit(FatalIOError);
1683 }
1684 }
1685 // Collect per face information
1686 label nExtrudeFaces = 0;
1687 forAll(meshZoneID, i)
1688 {
1689 nExtrudeFaces += faceZones[meshZoneID[i]].size();
1690 }
1691 extrudeMeshFaces.setSize(nExtrudeFaces);
1692 zoneFaces.setSize(nExtrudeFaces);
1693 zoneID.setSize(nExtrudeFaces);
1694 zoneFlipMap.setSize(nExtrudeFaces);
1695 nExtrudeFaces = 0;
1696 forAll(meshZoneID, i)
1697 {
1698 const faceZone& fz = faceZones[meshZoneID[i]];
1699 const primitiveFacePatch& fzp = fz();
1700 forAll(fz, j)
1701 {
1702 extrudeMeshFaces[nExtrudeFaces] = fz[j];
1703 zoneFaces[nExtrudeFaces] = fzp[j];
1704 zoneID[nExtrudeFaces] = i;
1705 zoneFlipMap[nExtrudeFaces] = fz.flipMap()[j];
1706 nExtrudeFaces++;
1707
1708 if (mesh.isInternalFace(fz[j]))
1709 {
1710 isInternal[i] = true;
1711 }
1712 }
1713 }
1714
1715 // Shadow zone
1716 // ~~~~~~~~~~~
1717
1718 if (zoneShadowNames.size())
1719 {
1720 zoneShadowIDs.setSize(zoneShadowNames.size());
1721 forAll(zoneShadowNames, i)
1722 {
1723 zoneShadowIDs[i] = faceZones.findZoneID(zoneShadowNames[i]);
1724 if (zoneShadowIDs[i] == -1)
1725 {
1727 << "Cannot find zone " << zoneShadowNames[i] << endl
1728 << "Valid zones are " << faceZones.names()
1729 << exit(FatalIOError);
1730 }
1731 }
1732
1733 label nShadowFaces = 0;
1734 forAll(zoneShadowIDs, i)
1735 {
1736 nShadowFaces += faceZones[zoneShadowIDs[i]].size();
1737 }
1738
1739 extrudeMeshShadowFaces.setSize(nShadowFaces);
1740 zoneShadowFlipMap.setSize(nShadowFaces);
1741 zoneShadowID.setSize(nShadowFaces);
1742
1743 nShadowFaces = 0;
1744 forAll(zoneShadowIDs, i)
1745 {
1746 const faceZone& fz = faceZones[zoneShadowIDs[i]];
1747 forAll(fz, j)
1748 {
1749 extrudeMeshShadowFaces[nShadowFaces] = fz[j];
1750 zoneShadowFlipMap[nShadowFaces] = fz.flipMap()[j];
1751 zoneShadowID[nShadowFaces] = i;
1752 nShadowFaces++;
1753 }
1754 }
1755 }
1756 }
1757 else
1758 {
1759 meshZoneID.setSize(zoneNames.size(), -1);
1760 // Load faceSets
1761 PtrList<faceSet> zones(zoneNames.size());
1762 forAll(zoneNames, i)
1763 {
1764 Info<< "Loading faceSet " << zoneNames[i] << endl;
1765 zones.set(i, new faceSet(mesh, zoneNames[i]));
1766 }
1767
1768
1769 // Collect per face information
1770 label nExtrudeFaces = 0;
1771 forAll(zones, i)
1772 {
1773 nExtrudeFaces += zones[i].size();
1774 }
1775 extrudeMeshFaces.setSize(nExtrudeFaces);
1776 zoneFaces.setSize(nExtrudeFaces);
1777 zoneID.setSize(nExtrudeFaces);
1778 zoneFlipMap.setSize(nExtrudeFaces);
1779
1780 nExtrudeFaces = 0;
1781 forAll(zones, i)
1782 {
1783 const faceSet& fz = zones[i];
1784 for (const label facei : fz)
1785 {
1786 if (mesh.isInternalFace(facei))
1787 {
1789 << "faceSet " << fz.name()
1790 << "contains internal faces."
1791 << " This is not permitted."
1792 << exit(FatalIOError);
1793 }
1794 extrudeMeshFaces[nExtrudeFaces] = facei;
1795 zoneFaces[nExtrudeFaces] = mesh.faces()[facei];
1796 zoneID[nExtrudeFaces] = i;
1797 zoneFlipMap[nExtrudeFaces] = false;
1798 nExtrudeFaces++;
1799
1800 if (mesh.isInternalFace(facei))
1801 {
1802 isInternal[i] = true;
1803 }
1804 }
1805 }
1806
1807
1808 // Shadow zone
1809 // ~~~~~~~~~~~
1810
1811 PtrList<faceSet> shadowZones(zoneShadowNames.size());
1812 if (zoneShadowNames.size())
1813 {
1814 zoneShadowIDs.setSize(zoneShadowNames.size(), -1);
1815 forAll(zoneShadowNames, i)
1816 {
1817 shadowZones.set(i, new faceSet(mesh, zoneShadowNames[i]));
1818 }
1819
1820 label nShadowFaces = 0;
1821 for (const faceSet& fz : shadowZones)
1822 {
1823 nShadowFaces += fz.size();
1824 }
1825
1826 if (nExtrudeFaces != nShadowFaces)
1827 {
1829 << "Extruded faces " << nExtrudeFaces << endl
1830 << "is different from shadow faces. " << nShadowFaces
1831 << "This is not permitted " << endl
1832 << exit(FatalIOError);
1833 }
1834
1835 extrudeMeshShadowFaces.setSize(nShadowFaces);
1836 zoneShadowFlipMap.setSize(nShadowFaces);
1837 zoneShadowID.setSize(nShadowFaces);
1838
1839 nShadowFaces = 0;
1840 forAll(shadowZones, i)
1841 {
1842 const faceSet& fz = shadowZones[i];
1843 for (const label facei : fz)
1844 {
1845 if (mesh.isInternalFace(facei))
1846 {
1848 << "faceSet " << fz.name()
1849 << "contains internal faces."
1850 << " This is not permitted."
1851 << exit(FatalIOError);
1852 }
1853 extrudeMeshShadowFaces[nShadowFaces] = facei;
1854 zoneShadowFlipMap[nShadowFaces] = false;
1855 zoneShadowID[nShadowFaces] = i;
1856 nShadowFaces++;
1857 }
1858 }
1859 }
1860 }
1861 const primitiveFacePatch extrudePatch(std::move(zoneFaces), mesh.points());
1862
1863
1865
1866 // Check zone either all internal or all external faces
1867 checkZoneInside(mesh, zoneNames, zoneID, extrudeMeshFaces, isInternal);
1868
1869
1870 const pointField& extrudePoints = extrudePatch.localPoints();
1871 const faceList& extrudeFaces = extrudePatch.localFaces();
1872 const labelListList& edgeFaces = extrudePatch.edgeFaces();
1873
1874
1875 Info<< "extrudePatch :"
1876 << " faces:" << extrudePatch.size()
1877 << " points:" << extrudePatch.nPoints()
1878 << " edges:" << extrudePatch.nEdges()
1879 << nl
1880 << endl;
1881
1882
1883 // Determine per-extrude-edge info
1884 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1885
1886 // Corresponding mesh edges
1887 const labelList extrudeMeshEdges
1888 (
1889 extrudePatch.meshEdges
1890 (
1891 mesh.edges(),
1892 mesh.pointEdges()
1893 )
1894 );
1895
1896 const globalIndex globalExtrudeFaces(extrudePatch.size());
1897
1898 // Global pp faces per pp edge.
1899 labelListList extrudeEdgeGlobalFaces
1900 (
1901 globalEdgeFaces
1902 (
1903 mesh,
1904 globalExtrudeFaces,
1905 extrudePatch,
1906 extrudeMeshEdges
1907 )
1908 );
1909 List<Map<label>> compactMap;
1910 const mapDistribute extrudeEdgeFacesMap
1911 (
1912 globalExtrudeFaces,
1913 extrudeEdgeGlobalFaces,
1914 compactMap
1915 );
1916
1917
1918 // Determine min and max zone per edge
1919 labelList edgeMinZoneID;
1920 labelList edgeMaxZoneID;
1921 calcEdgeMinMaxZone
1922 (
1923 mesh,
1924 extrudePatch,
1925 extrudeMeshEdges,
1926 zoneID,
1927 extrudeEdgeFacesMap,
1928 extrudeEdgeGlobalFaces,
1929
1930 edgeMinZoneID,
1931 edgeMaxZoneID
1932 );
1933
1934 if (!addSidePatches)
1935 {
1936 // Disabling inter-zone multiple boundary faces by setting 'left' and
1937 // 'right' zone to be the same
1938 edgeMaxZoneID = edgeMinZoneID;
1939 }
1940
1941
1942 DynamicList<polyPatch*> regionPatches(patches.size());
1943 // Copy all non-local patches since these are used on boundary edges of
1944 // the extrusion
1945 forAll(patches, patchi)
1946 {
1947 if (!isA<processorPolyPatch>(patches[patchi]))
1948 {
1949 label newPatchi = regionPatches.size();
1950 regionPatches.append
1951 (
1952 patches[patchi].clone
1953 (
1954 patches,
1955 newPatchi,
1956 0, // size
1957 0 // start
1958 ).ptr()
1959 );
1960 }
1961 }
1962
1963
1964 // Add interface patches
1965 // ~~~~~~~~~~~~~~~~~~~~~
1966
1967 // From zone to interface patch (region side)
1968 labelList interRegionTopPatch;
1969 labelList interRegionBottomPatch;
1970
1971 addCouplingPatches
1972 (
1973 mesh,
1974 regionName,
1975 shellRegionName,
1976 zoneNames,
1977 zoneShadowNames,
1978 isInternal,
1979 meshZoneID,
1980
1981 regionPatches,
1982 interRegionTopPatch,
1983 interRegionBottomPatch
1984 );
1985
1986
1987 // From zone to interface patch (mesh side)
1988 labelList interMeshTopPatch;
1989 labelList interMeshBottomPatch;
1990
1991 if (adaptMesh)
1992 {
1993 // Add coupling patches to mesh
1994
1995 // 1. Clone existing global patches
1996 DynamicList<polyPatch*> newPatches(patches.size());
1997 forAll(patches, patchi)
1998 {
1999 if (!isA<processorPolyPatch>(patches[patchi]))
2000 {
2001 autoPtr<polyPatch> clonedPatch(patches[patchi].clone(patches));
2002 clonedPatch->index() = newPatches.size();
2003 newPatches.append(clonedPatch.ptr());
2004 }
2005 }
2006
2007 // 2. Add new patches
2008 addCouplingPatches
2009 (
2010 mesh,
2011 regionName,
2012 shellRegionName,
2013 zoneNames,
2014 zoneShadowNames,
2015 isInternal,
2016 meshZoneID,
2017
2018 newPatches,
2019 interMeshTopPatch,
2020 interMeshBottomPatch
2021 );
2022
2023 // 3. Clone processor patches
2024 forAll(patches, patchi)
2025 {
2026 if (isA<processorPolyPatch>(patches[patchi]))
2027 {
2028 autoPtr<polyPatch> clonedPatch(patches[patchi].clone(patches));
2029 clonedPatch->index() = newPatches.size();
2030 newPatches.append(clonedPatch.ptr());
2031 }
2032 }
2033
2034
2035 #ifdef FULLDEBUG
2036 Pout<< "*** adaptMesh : addFvPAtches:" << nl;
2037 printPatches(Pout, newPatches);
2038 Pout<< "*** end of adaptMesh : addFvPAtches:" << endl;
2039 #endif
2040
2041
2042 // Add to mesh
2043 mesh.clearOut();
2044 mesh.removeFvBoundary();
2045 mesh.addFvPatches(newPatches, true);
2046
2048 }
2049
2050
2051 // Patch per extruded face
2052 labelList extrudeTopPatchID(extrudePatch.size());
2053 labelList extrudeBottomPatchID(extrudePatch.size());
2054
2055 forAll(zoneID, facei)
2056 {
2057 extrudeTopPatchID[facei] = interRegionTopPatch[zoneID[facei]];
2058 extrudeBottomPatchID[facei] = interRegionBottomPatch[zoneID[facei]];
2059 }
2060
2061
2062
2063 // Count how many patches on special edges of extrudePatch are necessary
2064 // - zoneXXX_sides
2065 // - zoneXXX_zoneYYY
2066 labelList zoneSidePatch(zoneNames.size(), Zero);
2067 // Patch to use for minZone
2068 labelList zoneZonePatch_min(zoneNames.size()*zoneNames.size(), Zero);
2069 // Patch to use for maxZone
2070 labelList zoneZonePatch_max(zoneNames.size()*zoneNames.size(), Zero);
2071
2072 countExtrudePatches
2073 (
2074 mesh,
2075 zoneNames.size(),
2076
2077 extrudePatch, // patch
2078 extrudeMeshFaces, // mesh face per patch face
2079 extrudeMeshEdges, // mesh edge per patch edge
2080
2081 extrudeEdgeGlobalFaces, // global indexing per patch edge
2082 edgeMinZoneID, // minZone per patch edge
2083 edgeMaxZoneID, // maxZone per patch edge
2084
2085 zoneSidePatch, // per zone-side num edges that extrude into it
2086 zoneZonePatch_min // per zone-zone num edges that extrude into it
2087 );
2088
2089 // Now we'll have:
2090 // zoneSidePatch[zoneA] : number of faces needed on the side of zoneA
2091 // zoneZonePatch_min[zoneA,zoneB] : number of faces needed inbetween A,B
2092
2093
2094 // Add the zone-side patches.
2095 addZoneSidePatches
2096 (
2097 mesh,
2098 zoneNames,
2099 (oneD ? oneDPatchType : word::null),
2100
2101 regionPatches,
2102 zoneSidePatch
2103 );
2104
2105
2106 // Add the patches inbetween zones
2107 addInterZonePatches
2108 (
2109 mesh,
2110 zoneNames,
2111 oneD,
2112
2113 zoneZonePatch_min,
2114 zoneZonePatch_max,
2115 regionPatches
2116 );
2117
2118
2119 // Sets sidePatchID[edgeI] to interprocessor patch. Adds any
2120 // interprocessor or cyclic patches if necessary.
2121 labelList sidePatchID;
2122 addCoupledPatches
2123 (
2124 mesh,
2125 extrudePatch,
2126 extrudeMeshFaces,
2127 extrudeMeshEdges,
2128 extrudeEdgeFacesMap,
2129 extrudeEdgeGlobalFaces,
2130
2131 sidePatchID,
2132 regionPatches
2133 );
2134
2135
2136// // Add all the newPatches to the mesh and fields
2137// // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2138// {
2139// forAll(newPatches, patchi)
2140// {
2141// Pout<< "Adding patch " << patchi
2142// << " name:" << newPatches[patchi]->name()
2143// << endl;
2144// }
2145// //label nOldPatches = mesh.boundary().size();
2146// mesh.clearOut();
2147// mesh.removeFvBoundary();
2148// mesh.addFvPatches(newPatches, true);
2149// //// Add calculated fvPatchFields for the added patches
2150// //for
2151// //(
2152// // label patchi = nOldPatches;
2153// // patchi < mesh.boundary().size();
2154// // patchi++
2155// //)
2156// //{
2157// // Pout<< "ADDing calculated to patch " << patchi
2158// // << endl;
2159// // addCalculatedPatchFields(mesh);
2160// //}
2161// //Pout<< "** Added " << mesh.boundary().size()-nOldPatches
2162// // << " patches." << endl;
2163// }
2164
2165
2166 // Set patches to use for edges to be extruded into boundary faces
2167 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2168 // In order of edgeFaces: per edge, per originating face the
2169 // patch to use for the side face (from the extruded edge).
2170 // If empty size create an internal face.
2171 labelListList extrudeEdgePatches(extrudePatch.nEdges());
2172
2173 // Is edge a non-manifold edge
2174 bitSet nonManifoldEdge(extrudePatch.nEdges());
2175
2176 // Note: logic has to be same as in countExtrudePatches.
2177 forAll(edgeFaces, edgeI)
2178 {
2179 const labelList& eFaces = edgeFaces[edgeI];
2180
2181 labelList& ePatches = extrudeEdgePatches[edgeI];
2182
2183 if (oneD)
2184 {
2185 ePatches.setSize(eFaces.size());
2186 forAll(eFaces, i)
2187 {
2188 ePatches[i] = zoneSidePatch[zoneID[eFaces[i]]];
2189 }
2190
2191 if (oneDNonManifoldEdges)
2192 {
2193 //- Set nonManifoldEdge[edgeI] for non-manifold edges only
2194 // The other option is to have non-manifold edges everywhere
2195 // and generate space overlapping columns of cells.
2196 if (eFaces.size() != 2)
2197 {
2198 nonManifoldEdge.set(edgeI);
2199 }
2200 }
2201 else
2202 {
2203 nonManifoldEdge.set(edgeI);
2204 }
2205 }
2206 else if (eFaces.size() == 2)
2207 {
2208 label zone0 = zoneID[eFaces[0]];
2209 label zone1 = zoneID[eFaces[1]];
2210
2211 if (addSidePatches && (zone0 != zone1)) // || (cos(angle) > blabla))
2212 {
2213 label minZone = Foam::min(zone0,zone1);
2214 label maxZone = Foam::max(zone0,zone1);
2215 label index = minZone*zoneNames.size()+maxZone;
2216
2217 ePatches.setSize(eFaces.size());
2218
2219 if (zone0 == minZone)
2220 {
2221 ePatches[0] = zoneZonePatch_min[index];
2222 ePatches[1] = zoneZonePatch_max[index];
2223 }
2224 else
2225 {
2226 ePatches[0] = zoneZonePatch_max[index];
2227 ePatches[1] = zoneZonePatch_min[index];
2228 }
2229
2230 nonManifoldEdge.set(edgeI);
2231 }
2232 }
2233 else if (sidePatchID[edgeI] != -1)
2234 {
2235 // Coupled extrusion
2236 ePatches.setSize(eFaces.size());
2237 forAll(eFaces, i)
2238 {
2239 ePatches[i] = sidePatchID[edgeI];
2240 }
2241 }
2242 else
2243 {
2244 label facei = findUncoveredPatchFace
2245 (
2246 mesh,
2247 labelUIndList(extrudeMeshFaces, eFaces),
2248 extrudeMeshEdges[edgeI]
2249 );
2250
2251 if (facei != -1)
2252 {
2253 label newPatchi = findPatchID
2254 (
2255 regionPatches,
2256 patches[patches.whichPatch(facei)].name()
2257 );
2258 ePatches.setSize(eFaces.size(), newPatchi);
2259 }
2260 else
2261 {
2262 ePatches.setSize(eFaces.size());
2263 forAll(eFaces, i)
2264 {
2265 ePatches[i] = zoneSidePatch[zoneID[eFaces[i]]];
2266 }
2267 }
2268 nonManifoldEdge.set(edgeI);
2269 }
2270 }
2271
2272
2273
2274 // Assign point regions
2275 // ~~~~~~~~~~~~~~~~~~~~
2276
2277 // Per face, per point the region number.
2278 faceList pointGlobalRegions;
2279 faceList pointLocalRegions;
2280 labelList localToGlobalRegion;
2281
2283 (
2284 mesh.globalData(),
2285 extrudePatch,
2286 nonManifoldEdge,
2287 false, // keep cyclic separated regions apart
2288
2289 pointGlobalRegions,
2290 pointLocalRegions,
2291 localToGlobalRegion
2292 );
2293
2294 // Per local region an originating point
2295 labelList localRegionPoints(localToGlobalRegion.size());
2296 forAll(pointLocalRegions, facei)
2297 {
2298 const face& f = extrudePatch.localFaces()[facei];
2299 const face& pRegions = pointLocalRegions[facei];
2300 forAll(pRegions, fp)
2301 {
2302 localRegionPoints[pRegions[fp]] = f[fp];
2303 }
2304 }
2305
2306 // Calculate region normals by reducing local region normals
2307 pointField localRegionNormals(localToGlobalRegion.size());
2308 {
2309 pointField localSum(localToGlobalRegion.size(), Zero);
2310
2311 forAll(pointLocalRegions, facei)
2312 {
2313 const face& pRegions = pointLocalRegions[facei];
2314 forAll(pRegions, fp)
2315 {
2316 label localRegionI = pRegions[fp];
2317 localSum[localRegionI] += extrudePatch.faceNormals()[facei];
2318 }
2319 }
2320
2321 Map<point> globalSum(2*localToGlobalRegion.size());
2322
2323 forAll(localSum, localRegionI)
2324 {
2325 label globalRegionI = localToGlobalRegion[localRegionI];
2326 globalSum.insert(globalRegionI, localSum[localRegionI]);
2327 }
2328
2329 // Reduce
2331
2332 forAll(localToGlobalRegion, localRegionI)
2333 {
2334 label globalRegionI = localToGlobalRegion[localRegionI];
2335 localRegionNormals[localRegionI] = globalSum[globalRegionI];
2336 }
2337 localRegionNormals /= mag(localRegionNormals);
2338 }
2339
2340
2341 // For debugging: dump hedgehog plot of normals
2342 if (false)
2343 {
2344 OFstream str(runTime.path()/"localRegionNormals.obj");
2345 label vertI = 0;
2346
2347 scalar thickness = model().sumThickness(1);
2348
2349 forAll(pointLocalRegions, facei)
2350 {
2351 const face& f = extrudeFaces[facei];
2352
2353 forAll(f, fp)
2354 {
2355 label region = pointLocalRegions[facei][fp];
2356 const point& pt = extrudePoints[f[fp]];
2357
2358 meshTools::writeOBJ(str, pt);
2359 vertI++;
2361 (
2362 str,
2363 pt+thickness*localRegionNormals[region]
2364 );
2365 vertI++;
2366 str << "l " << vertI-1 << ' ' << vertI << nl;
2367 }
2368 }
2369 }
2370
2371
2372 // Use model to create displacements of first layer
2373 vectorField firstDisp(localRegionNormals.size());
2374 forAll(firstDisp, regionI)
2375 {
2376 //const point& regionPt = regionCentres[regionI];
2377 const point& regionPt = extrudePatch.points()
2378 [
2379 extrudePatch.meshPoints()
2380 [
2381 localRegionPoints[regionI]
2382 ]
2383 ];
2384 const vector& n = localRegionNormals[regionI];
2385 firstDisp[regionI] = model()(regionPt, n, 1) - regionPt;
2386 }
2387
2388
2389 // Create a new mesh
2390 // ~~~~~~~~~~~~~~~~~
2391
2392 createShellMesh extruder
2393 (
2394 extrudePatch,
2395 pointLocalRegions,
2396 localRegionPoints
2397 );
2398
2399
2400 autoPtr<mapPolyMesh> shellMap;
2401 fvMesh regionMesh
2402 (
2403 IOobject
2404 (
2405 shellRegionName,
2406 meshInstance,
2407 runTime,
2411 ),
2412 Foam::zero{},
2413 false // syncPar
2414 );
2415
2416 // Add the new patches
2417 forAll(regionPatches, patchi)
2418 {
2419 polyPatch* ppPtr = regionPatches[patchi];
2420 regionPatches[patchi] = ppPtr->clone(regionMesh.boundaryMesh()).ptr();
2421 delete ppPtr;
2422 }
2423
2424 #ifdef FULLDEBUG
2425 Pout<< "*** regionPatches : regionPatches:" << nl;
2426 printPatches(Pout, regionPatches);
2427 Pout<< "*** end of regionPatches : regionPatches:" << endl;
2428 #endif
2429
2430
2431 regionMesh.clearOut();
2432 regionMesh.removeFvBoundary();
2433 regionMesh.addFvPatches(regionPatches, true);
2434
2435 {
2436 polyTopoChange meshMod(regionPatches.size());
2437
2438 extruder.setRefinement
2439 (
2440 firstDisp, // first displacement
2441 model().expansionRatio(),
2442 model().nLayers(), // nLayers
2443 extrudeTopPatchID,
2444 extrudeBottomPatchID,
2445 extrudeEdgePatches,
2446 meshMod
2447 );
2448
2449 // Enforce actual point positions according to extrudeModel (model)
2450 // (extruder.setRefinement only does fixed expansionRatio)
2451 // The regionPoints and nLayers are looped in the same way as in
2452 // createShellMesh
2453 DynamicList<point>& newPoints = const_cast<DynamicList<point>&>
2454 (
2455 meshMod.points()
2456 );
2457 label meshPointi = extrudePatch.localPoints().size();
2458 forAll(localRegionPoints, regionI)
2459 {
2460 label pointi = localRegionPoints[regionI];
2461 point pt = extrudePatch.localPoints()[pointi];
2462 const vector& n = localRegionNormals[regionI];
2463
2464 for (label layerI = 1; layerI <= model().nLayers(); layerI++)
2465 {
2466 newPoints[meshPointi++] = model()(pt, n, layerI);
2467 }
2468 }
2469
2470 shellMap = meshMod.changeMesh
2471 (
2472 regionMesh, // mesh to change
2473 false // inflate
2474 );
2475 }
2476
2477 // Necessary?
2478 regionMesh.setInstance(meshInstance);
2479
2480
2481 // Update numbering on extruder.
2482 extruder.updateMesh(shellMap());
2483
2484
2485 // Calculate offsets from shell mesh back to original mesh
2486 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2487
2488 List<pointField> topOffsets(zoneNames.size());
2489 List<boundBox> topBbs(zoneNames.size());
2490
2491 List<pointField> bottomOffsets(zoneNames.size());
2492 List<boundBox> bottomBbs(zoneNames.size());
2493
2494 forAll(regionMesh.boundaryMesh(), patchi)
2495 {
2496 const polyPatch& pp = regionMesh.boundaryMesh()[patchi];
2497
2499 {
2500 if (interRegionTopPatch.found(patchi))
2501 {
2502 label zoneI = interRegionTopPatch.find(patchi);
2503 topOffsets[zoneI] = calcOffset(extrudePatch, extruder, pp);
2504 topBbs[zoneI] = boundBox(pp.points(), pp.meshPoints(), true);
2505 }
2506 else if (interRegionBottomPatch.found(patchi))
2507 {
2508 label zoneI = interRegionBottomPatch.find(patchi);
2509 bottomOffsets[zoneI] = calcOffset(extrudePatch, extruder, pp);
2510 bottomBbs[zoneI] = boundBox(pp.points(), pp.meshPoints(), true);
2511 }
2512 }
2513 }
2514
2515
2516 // Change top and bottom boundary conditions on regionMesh
2517 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2518
2519 {
2520 // Correct top patches for offset
2521 setCouplingInfo
2522 (
2523 regionMesh,
2524 interRegionTopPatch,
2525 regionName, // name of main mesh
2526 sampleMode, // sampleMode
2527 topOffsets,
2528 topBbs
2529 );
2530
2531 // Correct bottom patches for offset
2532 setCouplingInfo
2533 (
2534 regionMesh,
2535 interRegionBottomPatch,
2536 regionName,
2537 sampleMode, // sampleMode
2538 bottomOffsets,
2539 bottomBbs
2540 );
2541
2542 // Remove any unused patches
2543 deleteEmptyPatches(regionMesh);
2544 }
2545
2546 // Change top and bottom boundary conditions on main mesh
2547 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2548
2549 if (adaptMesh)
2550 {
2551 // Correct top patches for offset
2552 setCouplingInfo
2553 (
2554 mesh,
2555 interMeshTopPatch,
2556 shellRegionName, // name of shell mesh
2557 sampleMode, // sampleMode
2558 -topOffsets,
2559 topBbs
2560 );
2561
2562 // Correct bottom patches for offset
2563 setCouplingInfo
2564 (
2565 mesh,
2566 interMeshBottomPatch,
2567 shellRegionName,
2568 sampleMode,
2569 -bottomOffsets,
2570 bottomBbs
2571 );
2572 }
2573
2574
2575
2576 // Write addressing from region mesh back to originating patch
2577 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2578
2579 labelIOList cellToPatchFaceAddressing
2580 (
2581 IOobject
2582 (
2583 "cellToPatchFaceAddressing",
2584 regionMesh.facesInstance(),
2586 regionMesh,
2590 ),
2591 extruder.cellToFaceMap()
2592 );
2593 cellToPatchFaceAddressing.note() = "cell to patch face addressing";
2594
2595 labelIOList faceToPatchFaceAddressing
2596 (
2597 IOobject
2598 (
2599 "faceToPatchFaceAddressing",
2600 regionMesh.facesInstance(),
2602 regionMesh,
2606 ),
2607 extruder.faceToFaceMap()
2608 );
2609 faceToPatchFaceAddressing.note() =
2610 "front/back face + turning index to patch face addressing";
2611
2612 labelIOList faceToPatchEdgeAddressing
2613 (
2614 IOobject
2615 (
2616 "faceToPatchEdgeAddressing",
2617 regionMesh.facesInstance(),
2619 regionMesh,
2623 ),
2624 extruder.faceToEdgeMap()
2625 );
2626 faceToPatchEdgeAddressing.note() =
2627 "side face to patch edge addressing";
2628
2629 labelIOList pointToPatchPointAddressing
2630 (
2631 IOobject
2632 (
2633 "pointToPatchPointAddressing",
2634 regionMesh.facesInstance(),
2636 regionMesh,
2640 ),
2641 extruder.pointToPointMap()
2642 );
2643 pointToPatchPointAddressing.note() =
2644 "point to patch point addressing";
2645
2646
2647 Info<< "Writing mesh " << regionMesh.name()
2648 << " to " << regionMesh.facesInstance() << nl
2649 << endl;
2650
2651 bool ok =
2652 regionMesh.write()
2653 && cellToPatchFaceAddressing.write()
2654 && faceToPatchFaceAddressing.write()
2655 && faceToPatchEdgeAddressing.write()
2656 && pointToPatchPointAddressing.write();
2657
2658 if (!ok)
2659 {
2661 << "Failed writing mesh " << regionMesh.name()
2662 << " at location " << regionMesh.facesInstance()
2663 << exit(FatalError);
2664 }
2665 topoSet::removeFiles(regionMesh);
2666 processorMeshes::removeFiles(regionMesh);
2667
2668
2669 // See if we need to extrude coordinates as well
2670 {
2671 autoPtr<pointIOField> patchFaceCentresPtr;
2672
2673 IOobject io
2674 (
2675 "patchFaceCentres",
2676 mesh.pointsInstance(),
2678 mesh,
2680 );
2681 if (io.typeHeaderOk<pointIOField>(true))
2682 {
2683 // Read patchFaceCentres and patchEdgeCentres
2684 Info<< "Reading patch face,edge centres : "
2685 << io.name() << " and patchEdgeCentres" << endl;
2686
2687 extrudeGeometricProperties
2688 (
2689 mesh,
2690 extrudePatch,
2691 extruder,
2692 regionMesh,
2693 model()
2694 );
2695 }
2696 }
2697
2698
2699
2700
2701 // Insert baffles into original mesh
2702 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2703
2704 autoPtr<mapPolyMesh> addBafflesMap;
2705
2706 if (adaptMesh)
2707 {
2708 polyTopoChange meshMod(mesh);
2709
2710 // Modify faces to be in bottom (= always coupled) patch
2711 forAll(extrudeMeshFaces, zoneFacei)
2712 {
2713 label meshFacei = extrudeMeshFaces[zoneFacei];
2714 label zoneI = zoneID[zoneFacei];
2715 bool flip = zoneFlipMap[zoneFacei];
2716 const face& f = mesh.faces()[meshFacei];
2717
2718 if (!flip)
2719 {
2720 meshMod.modifyFace
2721 (
2722 f, // modified face
2723 meshFacei, // label of face being modified
2724 mesh.faceOwner()[meshFacei],// owner
2725 -1, // neighbour
2726 false, // face flip
2727 interMeshBottomPatch[zoneI],// patch for face
2728 meshZoneID[zoneI], // zone for face
2729 flip // face flip in zone
2730 );
2731 }
2732 else if (mesh.isInternalFace(meshFacei))
2733 {
2734 meshMod.modifyFace
2735 (
2736 f.reverseFace(), // modified face
2737 meshFacei, // label of modified face
2738 mesh.faceNeighbour()[meshFacei],// owner
2739 -1, // neighbour
2740 true, // face flip
2741 interMeshBottomPatch[zoneI], // patch for face
2742 meshZoneID[zoneI], // zone for face
2743 !flip // face flip in zone
2744 );
2745 }
2746 }
2747
2748 if (zoneShadowNames.size() > 0) //if there is a top faceZone specified
2749 {
2750 forAll(extrudeMeshFaces, zoneFacei)
2751 {
2752 label meshFacei = extrudeMeshShadowFaces[zoneFacei];
2753 label zoneI = zoneShadowID[zoneFacei];
2754 bool flip = zoneShadowFlipMap[zoneFacei];
2755 const face& f = mesh.faces()[meshFacei];
2756
2757 if (!flip)
2758 {
2759 meshMod.modifyFace
2760 (
2761 f, // modified face
2762 meshFacei, // face being modified
2763 mesh.faceOwner()[meshFacei],// owner
2764 -1, // neighbour
2765 false, // face flip
2766 interMeshTopPatch[zoneI], // patch for face
2767 meshZoneID[zoneI], // zone for face
2768 flip // face flip in zone
2769 );
2770 }
2771 else if (mesh.isInternalFace(meshFacei))
2772 {
2773 meshMod.modifyFace
2774 (
2775 f.reverseFace(), // modified face
2776 meshFacei, // label modified face
2777 mesh.faceNeighbour()[meshFacei],// owner
2778 -1, // neighbour
2779 true, // face flip
2780 interMeshTopPatch[zoneI], // patch for face
2781 meshZoneID[zoneI], // zone for face
2782 !flip // face flip in zone
2783 );
2784 }
2785 }
2786 }
2787 else
2788 {
2789 // Add faces (using same points) to be in top patch
2790 forAll(extrudeMeshFaces, zoneFacei)
2791 {
2792 label meshFacei = extrudeMeshFaces[zoneFacei];
2793 label zoneI = zoneID[zoneFacei];
2794 bool flip = zoneFlipMap[zoneFacei];
2795 const face& f = mesh.faces()[meshFacei];
2796
2797 if (!flip)
2798 {
2799 if (mesh.isInternalFace(meshFacei))
2800 {
2801 meshMod.addFace
2802 (
2803 f.reverseFace(), // modified face
2804 mesh.faceNeighbour()[meshFacei],// owner
2805 -1, // neighbour
2806 -1, // master point
2807 -1, // master edge
2808 meshFacei, // master face
2809 true, // flip flux
2810 interMeshTopPatch[zoneI], // patch for face
2811 -1, // zone for face
2812 false //face flip in zone
2813 );
2814 }
2815 }
2816 else
2817 {
2818 meshMod.addFace
2819 (
2820 f, // face
2821 mesh.faceOwner()[meshFacei], // owner
2822 -1, // neighbour
2823 -1, // master point
2824 -1, // master edge
2825 meshFacei, // master face
2826 false, // flip flux
2827 interMeshTopPatch[zoneI], // patch for face
2828 -1, // zone for face
2829 false // zone flip
2830 );
2831 }
2832 }
2833 }
2834
2835 // Change the mesh. Change points directly (no inflation).
2836 addBafflesMap = meshMod.changeMesh(mesh, false);
2837
2838 // Update fields
2839 mesh.updateMesh(addBafflesMap());
2840
2841
2842//XXXXXX
2843// Update maps! e.g. faceToPatchFaceAddressing
2844//XXXXXX
2845
2846 // Move mesh (since morphing might not do this)
2847 if (addBafflesMap().hasMotionPoints())
2848 {
2849 mesh.movePoints(addBafflesMap().preMotionPoints());
2850 }
2851
2852 mesh.setInstance(meshInstance);
2853
2854 // Remove any unused patches
2855 deleteEmptyPatches(mesh);
2856
2857 Info<< "Writing mesh " << mesh.name()
2858 << " to " << mesh.facesInstance() << nl
2859 << endl;
2860
2861 if (!mesh.write())
2862 {
2864 << "Failed writing mesh " << mesh.name()
2865 << " at location " << mesh.facesInstance()
2866 << exit(FatalError);
2867 }
2870 }
2871
2872 Info << "End\n" << endl;
2873
2874 return 0;
2875}
2876
2877
2878// ************************************************************************* //
label n
Required Classes.
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
void append(const T &val)
Copy append an element to the end of this list.
SubField< vector > subField
Definition Field.H:183
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
label size() const noexcept
The number of elements in the list.
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
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 OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static tmp< pointField > edgeNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const labelList &patchEdges, const labelList &coupledEdges, const bitSet &flipMap=bitSet::null())
Return parallel consistent edge normals for patches using mesh points.
static void matchEdges(const PrimitivePatch< FaceList1, PointField1 > &p1, const PrimitivePatch< FaceList2, PointField2 > &p2, labelList &p1EdgeLabels, labelList &p2EdgeLabels, bitSet &sameOrientation)
Find corresponding edges on patches sharing the same points.
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & edgeFaces() const
Return edge-face addressing.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
static void mapCombineReduce(Container &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::mapReduce with an in-place cop.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void broadcast(Type &value, const int communicator=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
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.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
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 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
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
wordList names() const
A list of the zone names.
Definition ZoneMesh.C:463
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
static const Enum< transformType > transformTypeNames
Creates mesh by extruding a patch.
const labelList & faceToEdgeMap() const
From region side-face to patch edge. -1 for non-edge faces.
const labelList & pointToPointMap() const
From region point to patch point.
const labelList & faceToFaceMap() const
From region face to patch face. Contains turning index:
static void calcPointRegions(const globalMeshData &globalData, const primitiveFacePatch &patch, const bitSet &nonManifoldEdge, const bool syncNonCollocated, faceList &pointGlobalRegions, faceList &pointLocalRegions, labelList &localToGlobalRegion)
Helper: calculate point regions. The point region is the.
const labelList & cellToFaceMap() const
From region cell to patch face. Consecutively added so.
void updateMesh(const mapPolyMesh &)
Update any locally stored mesh information.
void setRefinement(const pointField &firstLayerThickness, const scalar expansionRatio, const label nLayers, const labelList &topPatchID, const labelList &bottomPatchID, const labelListList &extrudeEdgePatches, polyTopoChange &meshMod)
Play commands into polyTopoChange to create layer mesh.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
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
Top level extrusion model class.
static autoPtr< extrudeModel > New(const dictionary &dict)
Select null constructed.
label nLayers() const
Return the number of layers.
scalar sumThickness(const label layer) const
Helper: calculate cumulative relative thickness for layer.
A list of face labels.
Definition faceSet.H:50
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
const boolList & flipMap() const noexcept
Return face flip map.
Definition faceZone.H:389
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label toGlobal(const label proci, const label i) const
From local to global on proci.
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
sampleMode
Mesh items to sample.
static const Enum< sampleMode > sampleModeNames_
const Time & time() const noexcept
Return time registry.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition polyMeshIO.C:29
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition polyMesh.C:838
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
virtual autoPtr< polyPatch > clone(const labelList &faceCells) const
Construct and return a clone, setting faceCells.
Definition polyPatch.H:289
Direct mesh changes based on v1.3 polyTopoChange syntax.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const vectorField & faceCentres() const
const vectorField & cellCentres() const
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch ("procBoundary..") constructed from the pair of processor IDs...
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition topoSet.C:693
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
const polyBoundaryMesh & patches
bool coupled
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Required Classes.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
const word dictName("faMeshDefinition")
auto & name
const labelIOList & zoneIDs
Definition correctPhi.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
constexpr label labelMin
Definition label.H:54
UIndirectList< label > labelUIndList
UIndirectList of labels.
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition POSIX.C:775
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
vectorIOField pointIOField
pointIOField is a vectorIOField.
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
const pointField & pts
IOobject dictIO
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
List helper to append y unique elements onto the end of x.
Definition ListOps.H:721
Foam::surfaceFields.