Loading...
Searching...
No Matches
surfaceMeshExtract.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) 2017-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 surfaceMeshExtract
29
30Group
31 grpSurfaceUtilities
32
33Description
34 Extract patch or faceZone surfaces from a polyMesh.
35 Depending on output surface format triangulates faces.
36
37 Region numbers on faces not guaranteed to be the same as the patch indices.
38
39 Optionally only extracts named patches.
40
41 Optionally filters out points on faceZones, feature-edges and
42 featurePoints and generates pointPatches
43 for these - written to pointMesh/boundary.
44
45 If run in parallel, processor patches get filtered out by default and
46 the mesh is merged (based on topology).
47
48Usage
49 \b surfaceMeshExtract [OPTION] <surfacefile>
50
51 Options:
52 - \par -patches NAME | LIST
53 Specify single patch or multiple patches (name or regex) to extract
54
55 - \par -faceZones NAME | LIST
56 Specify single zone or multiple face zones (name or regex) to extract
57
58 - \par -exclude-patches NAME | LIST
59 Exclude single or multiple patches (name or regex) from extraction
60
61 - \par -excludeProcPatches
62 Exclude processor patches (default if parallel)
63
64 - \par -featureAngle <angle>
65 Extract feature edges/points and put into separate point-patches
66
67 - \par -extractZonePoints
68 Extract all face zone points and put into separate point-patches
69
70\*---------------------------------------------------------------------------*/
71
72#include "MeshedSurface.H"
74#include "argList.H"
75#include "Time.H"
76#include "polyMesh.H"
77#include "pointMesh.H"
78#include "emptyPolyPatch.H"
79#include "processorPolyPatch.H"
80#include "ListListOps.H"
81#include "stringListOps.H" // For stringListOps::findMatching()
83#include "globalMeshData.H"
84#include "globalIndex.H"
85#include "timeSelector.H"
86#include "meshPointPatch.H"
87#include "unitConversion.H"
88#include "dummyTransform.H"
89#include "syncTools.H"
90#include "processorPointPatch.H"
91#include "pointMeshTools.H"
92#include "OBJstream.H"
93
94using namespace Foam;
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98void writeOBJ
99(
100 const fileName& path,
101 const pointPatch& pp
102)
103{
105 const polyMesh& mesh = ppp.boundaryMesh().mesh().mesh();
106
107 // Count constraints
108 label maxConstraint = 0;
109 const auto& constraints = ppp.constraints();
110 forAll(constraints, i)
111 {
112 maxConstraint = Foam::max(maxConstraint, constraints[i].first());
113 }
114 reduce(maxConstraint, maxOp<label>());
115
116 const pointField localPoints(mesh.points(), ppp.meshPoints());
117
118 if (maxConstraint == 3)
119 {
121 (
122 path
123 / ppp.name()+"_fixedPoints.obj"
124 );
125 os.write(localPoints);
126 Info<< "Written pointPatch " << ppp.name() << " to " << os.name()
127 << endl;
128 }
129 else if (maxConstraint == 2)
130 {
132 (
133 path
134 / ppp.name()+"_slidingPoints.obj"
135 );
136 forAll(localPoints, i)
137 {
138 os.write(localPoints[i], constraints[i].second());
139 }
140 Info<< "Written pointPatch " << ppp.name() << " to " << os.name()
141 << " as coordinates and normals"
142 << endl;
143 }
144 else if (maxConstraint == 1)
145 {
147 (
148 path
149 / ppp.name()+"_surfacePoints.obj"
150 );
151 os.write(localPoints);
152 Info<< "Written pointPatch " << ppp.name() << " to " << os.name()
153 << " as coordinates"
154 << endl;
155 }
156}
157
158// Simple wrapper for polyBoundaryMesh::indices() with some additional logic
159// - prune emptyPolyPatch (always) and (maybe) processorPolyPatch
160labelList getSelectedPatches
161(
162 const polyBoundaryMesh& pbm,
163 const wordRes& allow,
164 const wordRes& deny,
165 const bool excludeProcPatches
166)
167{
168 labelList ids = pbm.indices(allow, deny);
169
170 // Remove undesirable patches
171 label count = 0;
172 for (const label patchi : ids)
173 {
174 const polyPatch& pp = pbm[patchi];
175
177 {
178 continue;
179 }
180 else if (excludeProcPatches && bool(isA<processorPolyPatch>(pp)))
181 {
182 break; // No processor patches for parallel output
183 }
184
185 ids[count] = patchi;
186 ++count;
187 }
188
189 ids.resize(count);
190 return ids;
191}
192
193
194label addMeshPointPatches
195(
196 const polyMesh& mesh,
197 const pointMesh& pMesh,
198
199 const uindirectPrimitivePatch& allBoundary,
200 const labelUList& faceToZone,
201 const surfZoneIdentifierList& surfZones,
202 const labelList& faceZoneToCompactZone,
203
204 const scalar edgeFeatureAngle,
205 const scalar pointFeatureAngle,
206 const bool verbose = true,
207 const bool writePoints = false
208)
209{
210 const auto& pointBm = pMesh.boundary();
211 const auto& fzs = mesh.faceZones();
212 const label nPointPatches = pointBm.size();
213 const pointField& points = mesh.points();
214
215
216 // Feature edge(points) internal to a zone
217 labelListList zoneToMeshPoints;
218 List<pointConstraintList> zoneToConstraints;
219
220 // Feature edge(points) in between zones
221 labelList twoZoneMeshPoints;
222 pointConstraintList twoZoneConstraints;
223
224 // Feature points on > 2 zones
225 labelList multiZoneMeshPoints;
226 pointConstraintList multiZoneConstraints;
227
229 (
230 mesh,
231 allBoundary,
232 // Per boundary face to zone
233 faceToZone,
234 // Number of zones
235 surfZones.size(),
236 edgeFeatureAngle,
237 //const scalar pointFeatureAngle, //not yet done
238
239 // Feature edge(points) internal to a zone
240 zoneToMeshPoints,
241 zoneToConstraints,
242
243 // Feature edge(points) in between zones
244 twoZoneMeshPoints,
245 twoZoneConstraints,
246
247 // Feature points on > 2 zones
248 multiZoneMeshPoints,
249 multiZoneConstraints
250 );
251
252
253 // Add per-zone patches
254 if (faceZoneToCompactZone.size())
255 {
256 // Calculate point normals consistent across whole patch
257 const pointField pointNormals
258 (
260 (
261 mesh,
262 allBoundary
263 )
264 );
265
266 forAll(faceZoneToCompactZone, zonei)
267 {
268 const label compacti = faceZoneToCompactZone[zonei];
269 if (compacti != -1)
270 {
271 const word patchName(surfZones[compacti].name());
272
273 if (pointBm.findPatchID(patchName) == -1)
274 {
275 // Extract the points originating from the faceZone. Can
276 // - re-call featurePointsEdges with 0 feature angle so
277 // all points go into the feature edges
278 // - mark using faceToZone the correct points
279 // - or assume the whole faceZone was extracted:
280 const uindirectPrimitivePatch fzPatch
281 (
283 (
284 mesh.faces(),
285 fzs[zonei].addressing()
286 ),
287 points
288 );
289 const auto& mp = fzPatch.meshPoints();
290
291 const vector nullVector(vector::uniform(0));
292
293 // Extract pointNormal (or 0) on all patch/connected points
294 vectorField meshPointNormals(mesh.nPoints(), nullVector);
295 for (const label pointi : mp)
296 {
297 const label allPointi =
298 allBoundary.meshPointMap()[pointi];
299 meshPointNormals[pointi] = pointNormals[allPointi];
300 }
302 (
303 mesh,
304 meshPointNormals,
306 nullVector
307 );
308
309 // Extract indices with non-zero pointNormal
310 DynamicList<label> meshPoints(mp.size());
311 forAll(meshPointNormals, pointi)
312 {
313 if (meshPointNormals[pointi] != nullVector)
314 {
315 meshPoints.append(pointi);
316 }
317 }
318
319 const_cast<pointBoundaryMesh&>(pointBm).push_back
320 (
322 (
323 patchName,
324 meshPoints,
325 vectorField(meshPointNormals, meshPoints),
326 pointBm.size(),
327 pointBm,
328 meshPointPatch::typeName
329 )
330 );
331
332 if (verbose)
333 {
334 const auto& ppp = pointBm.last();
335 Info<< "Added zone pointPatch " << ppp.name()
336 << " with "
337 << returnReduce(meshPoints.size(), sumOp<label>())
338 << " points" << endl;
339 }
340 if (writePoints)
341 {
342 writeOBJ(mesh.path(), pointBm.last());
343 }
344 }
345 }
346 }
347 }
348
349
350 // Add per-patch feature-edges
351 forAll(zoneToMeshPoints, zonei)
352 {
353 const label nPoints =
354 returnReduce(zoneToMeshPoints[zonei].size(), sumOp<label>());
355
356 const word patchName(surfZones[zonei].name() + "Edges");
357
358 if (nPoints && (pointBm.findPatchID(patchName) == -1))
359 {
360 const_cast<pointBoundaryMesh&>(pointBm).push_back
361 (
363 (
364 patchName,
365 zoneToMeshPoints[zonei],
366 zoneToConstraints[zonei],
367 pointBm.size(),
368 pointBm,
369 meshPointPatch::typeName
370 )
371 );
372
373 if (verbose)
374 {
375 const auto& ppp = pointBm.last();
376 Info<< "Added feature-edges pointPatch " << ppp.name()
377 << " with " << nPoints << " points" << endl;
378 }
379 if (writePoints)
380 {
381 writeOBJ(mesh.path(), pointBm.last());
382 }
383 }
384 }
385
386
387 // Add inter-patch points
388
389 const word allEdgePatchName("boundaryEdges");
390 const label nPatchEdgePoints =
391 returnReduce(twoZoneMeshPoints.size(), sumOp<label>());
392 if (nPatchEdgePoints && (pointBm.findPatchID(allEdgePatchName) == -1))
393 {
394 const_cast<pointBoundaryMesh&>(pointBm).push_back
395 (
397 (
398 allEdgePatchName,
399 twoZoneMeshPoints,
400 twoZoneConstraints,
401 pointBm.size(),
402 pointBm,
403 meshPointPatch::typeName
404 )
405 );
406
407 if (verbose)
408 {
409 const auto& ppp = pointBm.last();
410 Info<< "Added inter-patch pointPatch " << ppp.name()
411 << " with " << nPatchEdgePoints << " points" << endl;
412 }
413 if (writePoints)
414 {
415 writeOBJ(mesh.path(), pointBm.last());
416 }
417 }
418
419
420 const word allPointPatchName("boundaryPoints");
421 const label nMultiPoints =
422 returnReduce(multiZoneMeshPoints.size(), sumOp<label>());
423 if (nMultiPoints && (pointBm.findPatchID(allPointPatchName) == -1))
424 {
425 const_cast<pointBoundaryMesh&>(pointBm).push_back
426 (
428 (
429 allPointPatchName,
430 multiZoneMeshPoints,
431 multiZoneConstraints,
432 pointBm.size(),
433 pointBm,
434 meshPointPatch::typeName
435 )
436 );
437
438 if (verbose)
439 {
440 const auto& ppp = pointBm.last();
441 Info<< "Added multi-patch pointPatch " << ppp.name()
442 << " with " << nMultiPoints << " points" << endl;
443 }
444 if (writePoints)
445 {
446 writeOBJ(mesh.path(), pointBm.last());
447 }
448 }
449
450
451 // Shuffle into order
452 labelList oldToNew(pointBm.size());
453 label newPatchi = 0;
454 forAll(pointBm, patchi)
455 {
456 if (!isA<processorPointPatch>(pointBm[patchi]))
457 {
458 oldToNew[patchi] = newPatchi++;
459 }
460 }
461 forAll(pointBm, patchi)
462 {
463 if (isA<processorPointPatch>(pointBm[patchi]))
464 {
465 oldToNew[patchi] = newPatchi++;
466 }
467 }
468 const_cast<pointBoundaryMesh&>(pointBm).reorder(oldToNew, true);
469
470 return pointBm.size() - nPointPatches;
471}
472
473
474// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
475
476int main(int argc, char *argv[])
477{
479 (
480 "Extract patch or faceZone surfaces from a polyMesh."
481 );
483
484 // Less frequently used - reduce some clutter
485 argList::setAdvanced("decomposeParDict");
486
487 argList::addArgument("output", "The output surface file");
488
489 #include "addRegionOption.H"
491 (
492 "excludeProcPatches",
493 "Exclude processor patches"
494 );
496 (
497 "faceZones",
498 "wordRes",
499 "Specify single or multiple faceZones to extract\n"
500 "Eg, 'cells' or '( slice \"mfp-.*\" )'"
501 );
503 (
504 "patches",
505 "wordRes",
506 "Specify single patch or multiple patches to extract.\n"
507 "Eg, 'top' or '( front \".*back\" )'"
508 );
510 (
511 "exclude-patches",
512 "wordRes",
513 "Specify single patch or multiple patches to exclude from -patches."
514 " Eg, 'outlet' or '( inlet \".*Wall\" )'",
515 true // mark as an advanced option
516 );
517 argList::addOptionCompat("exclude-patches", {"excludePatches", 2306});
519 (
520 "featureAngle",
521 "angle",
522 "Auto-extract feature edges/points and put into separate point-patches"
523 );
525 (
526 "extractZonePoints",
527 "Extract point-patches for selected faceZones"
528 );
530 (
531 "writeOBJ",
532 "Write added pointPatch points to .obj files"
533 );
534
535 #include "setRootCase.H"
536 #include "createTime.H"
537
538 const auto userOutFileName = args.get<fileName>(1);
539
540 if (!userOutFileName.has_ext())
541 {
543 << "Missing extension on output name " << userOutFileName
544 << exit(FatalError);
545 }
546
547 Info<< "Extracting surface from boundaryMesh ..." << nl << nl;
548
549 const bool includeProcPatches =
550 (!UPstream::parRun() && !args.found("excludeProcPatches"));
551
552 if (includeProcPatches)
553 {
554 Info<< "Including all processor patches." << nl << endl;
555 }
556 else if (UPstream::parRun())
557 {
558 Info<< "Excluding all processor patches." << nl << endl;
559 }
560
561 wordRes includePatches, excludePatches;
562 if (args.readListIfPresent<wordRe>("patches", includePatches))
563 {
564 Info<< "Including patches " << flatOutput(includePatches)
565 << nl << endl;
566 }
567 if (args.readListIfPresent<wordRe>("exclude-patches", excludePatches))
568 {
569 Info<< "Excluding patches " << flatOutput(excludePatches)
570 << nl << endl;
571 }
572
573 // Non-mandatory
574 const wordRes selectedFaceZones(args.getList<wordRe>("faceZones", false));
575 if (selectedFaceZones.size())
576 {
577 Info<< "Including faceZones " << flatOutput(selectedFaceZones)
578 << nl << endl;
579 }
580
581 scalar featureAngle = 180.0;
582 const bool specifiedFeature = args.readIfPresent
583 (
584 "featureAngle",
585 featureAngle
586 );
587
588 const bool extractZonePoints = args.found("extractZonePoints");
589 const bool writeOBJ = args.found("writeOBJ");
590
591 Info<< "Reading mesh from time " << runTime.value() << endl;
592
593 #include "createNamedPolyMesh.H"
594 if (specifiedFeature)
595 {
596 Info<< "Detecting all sharp (>" << featureAngle
597 << " degrees) patch edges." << nl << endl;
598
599 if (extractZonePoints)
600 {
601 Info<< "Extracting all faceZone points as pointPatches."
602 << nl << endl;
603 }
604
605 //#include "createNamedPointMesh.H"
606 // Do not read constant/pointMesh - construct from polyMesh only
607 Info<< "Create pointMesh for time = "
608 << runTime.timeName() << Foam::nl << Foam::endl;
609 (void)pointMesh::New(mesh);
610 }
611
612
613 // User specified times
615
616 forAll(timeDirs, timeIndex)
617 {
618 runTime.setTime(timeDirs[timeIndex], timeIndex);
619 Info<< "Time [" << timeIndex << "] = " << runTime.timeName();
620
621 fileName outFileName;
622 if (timeDirs.size() == 1)
623 {
624 outFileName = userOutFileName;
625 }
626 else
627 {
630 {
631 Info<<" ... no mesh change." << nl;
632 continue;
633 }
634
635 // The filename based on the original, but with additional
636 // time information. The extension was previously checked that
637 // it exists
638 const auto dot = userOutFileName.rfind('.');
639
640 outFileName =
641 userOutFileName.substr(0, dot) + "_"
642 + Foam::name(runTime.value()) + "."
643 + userOutFileName.ext();
644 }
645
646 Info<< nl;
647
648 // Create local surface from:
649 // - explicitly named patches only (-patches (at your option)
650 // - all patches (default in sequential mode)
651 // - all non-processor patches (default in parallel mode)
652 // - all non-processor patches (sequential mode, -excludeProcPatches
653 // (at your option)
654
655 // Construct table of patches to include.
656 const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
657
658 const labelList patchIds =
659 (
660 getSelectedPatches
661 (
662 bMesh,
663 includePatches,
664 excludePatches,
665 // No processor patches? (parallel output or excluded)
666 (UPstream::parRun() || !includeProcPatches)
667 )
668 );
669
670 labelList faceZoneIds;
671
672 const faceZoneMesh& fzm = mesh.faceZones();
673
674 if (selectedFaceZones.size())
675 {
676 faceZoneIds = fzm.indices(selectedFaceZones);
677
678 Info<< "Additionally extracting faceZones "
679 << fzm.names(selectedFaceZones) << nl;
680 }
681
682
683 // From (name of) patch to compact 'zone' index
684 HashTable<label> compactZoneID(1024);
685 // Mesh face and compact zone indx
687 DynamicList<label> compactZones;
688 // Per compact 'zone' index the name and location
689 surfZoneIdentifierList surfZones;
690
691 // Per local patch to compact 'zone' index (or -1)
692 labelList patchToCompactZone(bMesh.size(), -1);
693 // Per local faceZone to compact 'zone' index (or -1)
694 labelList faceZoneToCompactZone(bMesh.size(), -1);
695 {
696 // Collect sizes. Hash on names to handle local-only patches (e.g.
697 // processor patches)
698 HashTable<label> patchSize(1024);
699 label nFaces = 0;
700 for (const label patchi : patchIds)
701 {
702 const polyPatch& pp = bMesh[patchi];
703 patchSize.insert(pp.name(), pp.size());
704 nFaces += pp.size();
705 }
706
707 HashTable<label> zoneSize(1024);
708 for (const label zonei : faceZoneIds)
709 {
710 const faceZone& pp = fzm[zonei];
711 zoneSize.insert(pp.name(), pp.size());
712 nFaces += pp.size();
713 }
714
715
718
719 if (Pstream::master())
720 {
721 // Allocate compact numbering for all patches/faceZones
722 forAllConstIters(patchSize, iter)
723 {
724 compactZoneID.insert(iter.key(), compactZoneID.size());
725 }
726
727 forAllConstIters(zoneSize, iter)
728 {
729 compactZoneID.insert(iter.key(), compactZoneID.size());
730 }
731 }
732 Pstream::broadcast(compactZoneID);
733
734
735 // Zones
736 surfZones.resize_nocopy(compactZoneID.size());
737 forAllConstIters(compactZoneID, iter)
738 {
739 surfZones[*iter] = surfZoneIdentifier(iter.key(), *iter);
740 Info<< "surfZone " << *iter
741 << " : " << surfZones[*iter].name()
742 << endl;
743 }
744
745
746 // Rework HashTable into labelList just for speed of conversion
747 forAllConstIters(compactZoneID, iter)
748 {
749 label patchi = bMesh.findPatchID(iter.key());
750 if (patchi != -1)
751 {
752 patchToCompactZone[patchi] = iter.val();
753 }
754 else
755 {
756 label zoneI = fzm.findZoneID(iter.key());
757 faceZoneToCompactZone[zoneI] = iter.val();
758 }
759 }
760
761
762 faceLabels.setCapacity(nFaces);
763 compactZones.setCapacity(nFaces);
764
765 // Collect faces on patches
766 for (const label patchi : patchIds)
767 {
768 const polyPatch& pp = bMesh[patchi];
769 forAll(pp, i)
770 {
771 faceLabels.append(pp.start()+i);
772 compactZones.append(patchToCompactZone[pp.index()]);
773 }
774 }
775 // Collect faces on faceZones
776 for (const label zonei : faceZoneIds)
777 {
778 const faceZone& pp = fzm[zonei];
779 forAll(pp, i)
780 {
781 faceLabels.append(pp[i]);
782 compactZones.append(faceZoneToCompactZone[pp.index()]);
783 }
784 }
785 }
786
787
788 // Addressing engine for all faces
789 const uindirectPrimitivePatch allBoundary
790 (
792 mesh.points()
793 );
794
795
796 // Find correspondence to master points
797 labelList pointToGlobal;
798 labelList uniqueMeshPoints;
799 autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
800 (
801 allBoundary.meshPoints(),
802 allBoundary.meshPointMap(),
803 pointToGlobal,
804 uniqueMeshPoints
805 );
806
807 // Gather all unique points on master
808 List<pointField> gatheredPoints(Pstream::nProcs());
809 gatheredPoints[Pstream::myProcNo()] = pointField
810 (
811 mesh.points(),
812 uniqueMeshPoints
813 );
814 Pstream::gatherList(gatheredPoints);
815
816 // Gather all faces
817 List<faceList> gatheredFaces(Pstream::nProcs());
818 gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
819 for (face& f : gatheredFaces[Pstream::myProcNo()])
820 {
821 inplaceRenumber(pointToGlobal, f);
822 }
823 Pstream::gatherList(gatheredFaces);
824
825 // Gather all ZoneIDs
826 List<labelList> gatheredZones(Pstream::nProcs());
827 gatheredZones[Pstream::myProcNo()] = compactZones;
828 Pstream::gatherList(gatheredZones);
829
830 // On master combine all points, faces, zones
831 if (Pstream::master())
832 {
834 (
835 gatheredPoints,
837 );
838 gatheredPoints.clear();
839
841 (
842 gatheredFaces,
844 );
845 gatheredFaces.clear();
846
848 (
849 gatheredZones,
851 );
852 gatheredZones.clear();
853
854
855 UnsortedMeshedSurface<face> unsortedFace
856 (
857 std::move(allPoints),
858 std::move(allFaces),
859 std::move(allZones),
860 surfZones
861 );
862
863
864 MeshedSurface<face> sortedFace(unsortedFace);
865
866 fileName globalCasePath
867 (
868 outFileName.isAbsolute()
869 ? outFileName
870 : (
871 runTime.processorCase()
872 ? runTime.globalPath()/outFileName
873 : runTime.path()/outFileName
874 )
875 );
876
877 Info<< "Writing merged surface to " << globalCasePath << endl;
878
879 sortedFace.write(globalCasePath);
880 }
881
882
883 if (specifiedFeature)
884 {
885 // Add edge patches
886 const auto& pMesh = pointMesh::New(mesh);
887
888 const label nAdded = addMeshPointPatches
889 (
890 mesh,
891 pMesh,
892
893 allBoundary, // all patches together
894 compactZones, // originating compactZone
895 surfZones, // per compactZone the index
896 (
897 extractZonePoints
898 ? faceZoneToCompactZone // per faceZone the compactZone
900 ),
901
902 featureAngle,
903 featureAngle,
904
905 true,
906 writeOBJ
907 );
908
909 if (nAdded)
910 {
911 pMesh.boundary().write();
912 }
913 }
914 }
915
916 Info<< "End\n" << endl;
917
918 return 0;
919}
920
921
922// ************************************************************************* //
Required Classes.
labelList faceLabels(nFaceLabels)
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.
void setCapacity(const label len)
Alter the size of the underlying storage.
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
static const List< label > & null() noexcept
Definition List.H:138
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
const Mesh & mesh() const noexcept
Reference to the mesh.
Definition MeshObject.H:257
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const bitSet &flipMap=bitSet::null())
Return parallel consistent point normals for patches using mesh points.
const Map< label > & meshPointMap() const
Mesh point map.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
static void mapCombineGather(Container &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::mapGather with an in-place cop.
static void gatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate.
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...
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
A surface geometry mesh, in which the surface zone information is conveyed by the 'zoneId' associated...
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
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
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching zone groups.
Definition ZoneMesh.C:562
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition argList.C:366
static void setAdvanced(const word &optName, bool advanced=true)
Set an existing option as being 'advanced' or normal.
Definition argList.C:419
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
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
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition argList.C:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling file names.
Definition fileName.H:75
static bool isAbsolute(const std::string &str)
Return true if filename starts with a '/' or '\' or (windows-only) with a filesystem-root.
Definition fileNameI.H:129
pointPatch with explicitly provided points instead of using the points of a polyPatch.
virtual const List< pointConstraint > & constraints() const
Return constraints.
virtual const labelList & meshPoints() const
Return mesh points.
Database for mesh data, solution data, solver performance and other reduced data.
Definition meshState.H:56
const word & name() const noexcept
The patch name.
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
const pointMesh & mesh() const noexcept
Return the mesh reference.
static void featurePointsEdges(const polyMesh &mesh, const uindirectPrimitivePatch &boundary, const labelUList &faceToZone, const label nZones, const scalar edgeFeatureAngle, labelListList &zoneToMeshPoints, List< pointConstraintList > &zoneToConstraints, labelList &twoZoneMeshPoints, pointConstraintList &twoZoneConstraints, labelList &multiZoneMeshPoints, pointConstraintList &multiZoneConstraints)
Analyse patch for feature edges, feature points. Handles points not being on a face of patch but coup...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
const pointBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition pointMesh.H:177
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
const pointBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition pointPatch.H:250
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
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition polyMesh.H:92
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Identifies a surface patch/zone by name and index, with optional geometric type.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
word ext() const
Return file name extension (part after last .).
Definition wordI.H:171
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
labelList patchIds
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
Dummy transform to be used with syncTools.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
label nPoints
return returnReduce(nRefine-oldNRefine, sumOp< label >())
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition ListListOps.C:62
const dimensionedScalar mp
Proton mass.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
List< surfZoneIdentifier > surfZoneIdentifierList
List of surfZoneIdentifier.
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< pointConstraint > pointConstraintList
List of pointConstraint.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
List< instant > instantList
List of instants.
Definition instantList.H:41
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
Field< vector > vectorField
Specialisation of Field<T> for vector.
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
ListType reorder(const labelUList &oldToNew, const ListType &input, const bool prune=false)
Reorder the elements of a list.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
labelList f(nPoints)
Foam::argList args(argc, argv)
#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
Operations on lists of strings.
Object access operator or list access operator (default is pass-through).
Definition UList.H:1120
Unit conversion functions.