Loading...
Searching...
No Matches
conformalVoronoiMeshIO.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "IOstreams.H"
31#include "OFstream.H"
32#include "pointMesh.H"
33#include "pointFields.H"
34#include "ListOps.H"
35#include "polyMeshFilter.H"
36#include "polyTopoChange.H"
37#include "PrintTable.H"
38#include "indexedVertexOps.H"
39#include "DelaunayMeshTools.H"
40#include "syncTools.H"
41#include "memInfo.H"
42#include "faceSet.H"
43#include "OBJstream.H"
44
45// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46
48(
49 const string& description
50) const
51{
52 timeCheck(time(), description, foamyHexMeshControls().timeChecks());
53}
54
55
57(
58 const Time& runTime,
59 const string& description,
60 const bool check
61)
62{
63 if (check)
64 {
65 Info<< nl << "--- [ cpuTime "
66 << runTime.elapsedCpuTime() << " s, "
67 << "delta " << runTime.cpuTimeIncrement()<< " s";
68
69 if (!description.empty())
70 {
71 Info<< ", " << description << " ";
72 }
73 else
74 {
75 Info<< " ";
76 }
77
78 Info<< "] --- " << endl;
79
80 memInfo m;
81
82 if (m.good())
83 {
84 PrintTable<word, label> memoryTable
85 (
86 "Memory Usage (kB): "
87 + description
88 );
89
90 memoryTable.add("mSize", m.size());
91 memoryTable.add("mPeak", m.peak());
92 memoryTable.add("mRss", m.rss());
93
95 memoryTable.print(Info, true, true);
97 }
98 }
99}
100
101
103{
105
106 // Per cell the Delaunay vertex
107 labelList cellToDelaunayVertex;
108 // Per patch, per face the Delaunay vertex
109 labelListList patchToDelaunayVertex;
110
111 labelList dualPatchStarts;
112
113 {
115 labelList boundaryPts;
116 faceList faces;
117 labelList owner;
118 labelList neighbour;
121 pointField cellCentres;
122 bitSet boundaryFacesToRemove;
123
124 calcDualMesh
125 (
126 points,
127 boundaryPts,
128 faces,
129 owner,
130 neighbour,
133 cellCentres,
134 cellToDelaunayVertex,
135 patchToDelaunayVertex,
136 boundaryFacesToRemove
137 );
138
139 Info<< nl << "Writing polyMesh to " << instance << endl;
140
141 writeMesh
142 (
144 instance,
145 points,
146 boundaryPts,
147 faces,
148 owner,
149 neighbour,
152 cellCentres,
153 boundaryFacesToRemove
154 );
155
156 dualPatchStarts.setSize(patchDicts.size());
157
158 forAll(dualPatchStarts, patchi)
159 {
160 dualPatchStarts[patchi] =
161 patchDicts[patchi].get<label>("startFace");
162 }
163 }
164
165 if (foamyHexMeshControls().writeCellShapeControlMesh())
166 {
167 cellShapeControls().shapeControlMesh().write();
168 }
169
170 if (foamyHexMeshControls().writeBackgroundMeshDecomposition())
171 {
172 Info<< nl << "Writing " << "backgroundMeshDecomposition" << endl;
173
174 // Have to explicitly update the mesh instance.
175 const_cast<fvMesh&>(decomposition_().mesh()).setInstance
176 (
177 time().timeName()
178 );
179
180 decomposition_().mesh().write();
181 }
182
183 if (foamyHexMeshControls().writeTetDualMesh())
184 {
185 label celli = 0;
186 for
187 (
188 Finite_cells_iterator cit = finite_cells_begin();
189 cit != finite_cells_end();
190 ++cit
191 )
192 {
193 if
194 (
195 !cit->hasFarPoint()
196 && !is_infinite(cit)
197 )
198 {
199 cit->cellIndex() = celli++;
200 }
201 }
202
203 Info<< nl << "Writing " << "tetDualMesh" << endl;
204
205 labelPairLookup vertexMap;
206 labelList cellMap;
207 autoPtr<polyMesh> tetMesh =
208 createMesh("tetDualMesh", vertexMap, cellMap);
209
210 tetMesh().write();
211
212// // Determine map from Delaunay vertex to Dual mesh
213// // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
214//
215// // From all Delaunay vertices to cell (positive index)
216// // or patch face (negative index)
217// labelList vertexToDualAddressing(number_of_vertices(), Zero);
218//
219// forAll(cellToDelaunayVertex, celli)
220// {
221// label vertI = cellToDelaunayVertex[celli];
222//
223// if (vertexToDualAddressing[vertI] != 0)
224// {
225// FatalErrorInFunction
226// << "Delaunay vertex " << vertI
227// << " from cell " << celli
228// << " is already mapped to "
229// << vertexToDualAddressing[vertI]
230// << exit(FatalError);
231// }
232// vertexToDualAddressing[vertI] = celli+1;
233// }
234//
235// forAll(patchToDelaunayVertex, patchi)
236// {
237// const labelList& patchVertices = patchToDelaunayVertex[patchi];
238//
239// forAll(patchVertices, i)
240// {
241// label vertI = patchVertices[i];
242//
243// if (vertexToDualAddressing[vertI] > 0)
244// {
245// FatalErrorInFunction
246// << "Delaunay vertex " << vertI
247// << " from patch " << patchi
248// << " local index " << i
249// << " is already mapped to cell "
250// << vertexToDualAddressing[vertI]-1
251// << exit(FatalError);
252// }
253//
254// // Vertex might be used by multiple faces. Which one to
255// // use? For now last one wins.
256// label dualFacei = dualPatchStarts[patchi]+i;
257// vertexToDualAddressing[vertI] = -dualFacei-1;
258// }
259// }
260//
261//
262// // Calculate tet mesh addressing
263// // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
264//
265// pointField points;
266// labelList boundaryPts(number_of_finite_cells(), -1);
267// // From tet point back to Delaunay vertex index
268// labelList pointToDelaunayVertex;
269// faceList faces;
270// labelList owner;
271// labelList neighbour;
272// wordList patchTypes;
273// wordList patchNames;
274// PtrList<dictionary> patchDicts;
275// pointField cellCentres;
276//
277// calcTetMesh
278// (
279// points,
280// pointToDelaunayVertex,
281// faces,
282// owner,
283// neighbour,
284// patchTypes,
285// patchNames,
286// patchDicts
287// );
288//
289//
290//
291// // Calculate map from tet points to dual mesh cells/patch faces
292// // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
293//
294// labelIOList pointDualAddressing
295// (
296// IOobject
297// (
298// "pointDualAddressing",
299// instance,
300// polyMesh::meshDir("tetDualMesh"),
301// runTime_,
302// IOobject::NO_READ,
303// IOobject::NO_WRITE,
304// IOobject::NO_REGISTER
305// ),
306// labelUIndList
307// (
308// vertexToDualAddressing,
309// pointToDelaunayVertex
310// )()
311// );
312//
313// label pointi = pointDualAddressing.find(-1);
314// if (pointi != -1)
315// {
316// WarningInFunction
317// << "Delaunay vertex " << pointi
318// << " does not have a corresponding dual cell." << endl;
319// }
320//
321// Info<< "Writing map from tetDualMesh points to Voronoi mesh to "
322// << pointDualAddressing.objectPath() << endl;
323// pointDualAddressing.write();
324//
325//
326//
327// // Write tet points corresponding to the Voronoi cell/face centre
328// // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
329// {
330// // Read Voronoi mesh
331// fvMesh mesh
332// (
333// IOobject
334// (
335// polyMesh::defaultRegion,
336// instance,
337// runTime_,
338// IOobject::MUST_READ
339// )
340// );
341// pointIOField dualPoints
342// (
343// IOobject
344// (
345// "dualPoints",
346// instance,
347// polyMesh::meshDir("tetDualMesh"),
348// runTime_,
349// IOobject::NO_READ,
350// IOobject::NO_WRITE,
351// IOobject::NO_REGISTER
352// ),
353// points
354// );
355//
356// forAll(pointDualAddressing, pointi)
357// {
358// label index = pointDualAddressing[pointi];
359//
360// if (index > 0)
361// {
362// label celli = index-1;
363// dualPoints[pointi] = mesh.cellCentres()[celli];
364// }
365// else if (index < 0)
366// {
367// label facei = -index-1;
368// if (facei >= mesh.nInternalFaces())
369// {
370// dualPoints[pointi] = mesh.faceCentres()[facei];
371// }
372// }
373// }
374//
375// Info<< "Writing tetDualMesh points mapped onto Voronoi mesh to "
376// << dualPoints.objectPath() << endl
377// << "Replace the polyMesh/points with these." << endl;
378// dualPoints.write();
379// }
380//
381//
382// Info<< nl << "Writing tetDualMesh to " << instance << endl;
383//
384// bitSet boundaryFacesToRemove;
385// writeMesh
386// (
387// "tetDualMesh",
388// instance,
389// points,
390// boundaryPts,
391// faces,
392// owner,
393// neighbour,
394// patchTypes,
395// patchNames,
396// patchDicts,
397// cellCentres,
398// boundaryFacesToRemove
399// );
400 }
401}
402
403
404Foam::autoPtr<Foam::fvMesh> Foam::conformalVoronoiMesh::createDummyMesh
405(
406 const IOobject& io,
407 const wordList& patchNames,
409) const
410{
412 fvMesh& mesh = meshPtr();
413
415
416 forAll(patches, patchi)
417 {
418 if
419 (
420 patchDicts.set(patchi)
421 && (
422 patchDicts[patchi].get<word>("type")
423 == processorPolyPatch::typeName
424 )
425 )
426 {
427 patches[patchi] = new processorPolyPatch
428 (
429 0, //patchSizes[p],
430 0, //patchStarts[p],
431 patchi,
433 patchDicts[patchi].get<label>("myProcNo"),
434 patchDicts[patchi].get<label>("neighbProcNo"),
436 );
437 }
438 else
439 {
440 patches[patchi] = polyPatch::New
441 (
442 patchDicts[patchi].get<word>("type"),
443 patchNames[patchi],
444 0, //patchSizes[p],
445 0, //patchStarts[p],
446 patchi,
448 ).ptr();
449 }
450 }
451
453
454 return meshPtr;
455}
456
457
458void Foam::conformalVoronoiMesh::checkProcessorPatchesMatch
459(
461) const
462{
463 // Check patch sizes
464 labelListList procPatchSizes
465 (
468 );
469
470 forAll(patchDicts, patchi)
471 {
472 if
473 (
474 patchDicts.set(patchi)
475 && (
476 patchDicts[patchi].get<word>("type")
477 == processorPolyPatch::typeName
478 )
479 )
480 {
481 const label procNeighb =
482 patchDicts[patchi].get<label>("neighbProcNo");
483
484 procPatchSizes[Pstream::myProcNo()][procNeighb]
485 = patchDicts[patchi].get<label>("nFaces");
486 }
487 }
488
489 Pstream::gatherList(procPatchSizes);
490
491 if (Pstream::master())
492 {
493 bool allMatch = true;
494
495 forAll(procPatchSizes, proci)
496 {
497 const labelList& patchSizes = procPatchSizes[proci];
498
499 forAll(patchSizes, patchi)
500 {
501 if (patchSizes[patchi] != procPatchSizes[patchi][proci])
502 {
503 allMatch = false;
504
505 Info<< indent << "Patches " << proci << " and " << patchi
506 << " have different sizes: " << patchSizes[patchi]
507 << " and " << procPatchSizes[patchi][proci] << endl;
508 }
509 }
510 }
511
512 if (allMatch)
513 {
514 Info<< indent << "All processor patches have matching numbers of "
515 << "faces" << endl;
516 }
517 }
518}
519
520
521void Foam::conformalVoronoiMesh::reorderPoints
522(
524 labelList& boundaryPts,
525 faceList& faces,
526 const label nInternalFaces
527) const
528{
529 Info<< incrIndent << indent << "Reordering points into internal/external"
530 << endl;
531
532 labelList oldToNew(points.size(), Zero);
533
534 // Find points that are internal
535 for (label fI = nInternalFaces; fI < faces.size(); ++fI)
536 {
537 const face& f = faces[fI];
538
539 forAll(f, fpI)
540 {
541 oldToNew[f[fpI]] = 1;
542 }
543 }
544
545 const label nInternalPoints = points.size() - sum(oldToNew);
546
547 label countInternal = 0;
548 label countExternal = nInternalPoints;
549
550 forAll(points, pI)
551 {
552 if (oldToNew[pI] == 0)
553 {
554 oldToNew[pI] = countInternal++;
555 }
556 else
557 {
558 oldToNew[pI] = countExternal++;
559 }
560 }
561
562 Info<< indent
563 << "Number of internal points: " << countInternal << nl
564 << indent << "Number of external points: " << countExternal
565 << decrIndent << endl;
566
567 inplaceReorder(oldToNew, points);
568 inplaceReorder(oldToNew, boundaryPts);
569
570 forAll(faces, fI)
571 {
572 face& f = faces[fI];
573
574 forAll(f, fpI)
575 {
576 f[fpI] = oldToNew[f[fpI]];
577 }
578 }
579}
580
581
582void Foam::conformalVoronoiMesh::reorderProcessorPatches
583(
584 const word& meshName,
585 const fileName& instance,
586 const pointField& points,
587 faceList& faces,
588 const wordList& patchNames,
590) const
591{
592 Info<< incrIndent << indent << "Reordering processor patches" << endl;
593
595 checkProcessorPatchesMatch(patchDicts);
596
597 // Create dummy mesh with correct proc boundaries to do sorting
598 autoPtr<fvMesh> sortMeshPtr
599 (
600 createDummyMesh
601 (
603 (
604 meshName,
605 instance,
606 runTime_,
610 ),
613 )
614 );
615 const fvMesh& sortMesh = sortMeshPtr();
616
617 // Rotation on new faces.
618 labelList rotation(faces.size(), Zero);
619 labelList faceMap(faces.size(), label(-1));
620
621 PstreamBuffers pBufs;
622
623 // Send ordering
624 forAll(sortMesh.boundaryMesh(), patchi)
625 {
626 const polyPatch& pp = sortMesh.boundaryMesh()[patchi];
627
629 {
631 (
632 pBufs,
634 (
636 (
637 faces,
638 patchDicts[patchi].get<label>("nFaces"),
639 patchDicts[patchi].get<label>("startFace")
640 ),
641 points
642 )
643 );
644 }
645 }
646
647 pBufs.finishedSends();
648
649 Info<< incrIndent << indent << "Face ordering initialised..." << endl;
650
651 // Receive and calculate ordering
652 bool anyChanged = false;
653
654 forAll(sortMesh.boundaryMesh(), patchi)
655 {
656 const polyPatch& pp = sortMesh.boundaryMesh()[patchi];
657
659 {
660 const label nPatchFaces =
661 patchDicts[patchi].get<label>("nFaces");
662
663 const label patchStartFace =
664 patchDicts[patchi].get<label>("startFace");
665
666 labelList patchFaceMap(nPatchFaces, label(-1));
667 labelList patchFaceRotation(nPatchFaces, Zero);
668
669 bool changed = refCast<const processorPolyPatch>(pp).order
670 (
671 pBufs,
673 (
675 (
676 faces,
677 nPatchFaces,
678 patchStartFace
679 ),
680 points
681 ),
682 patchFaceMap,
683 patchFaceRotation
684 );
685
686 if (changed)
687 {
688 // Merge patch face reordering into mesh face reordering table
689 forAll(patchFaceRotation, patchFacei)
690 {
691 rotation[patchFacei + patchStartFace]
692 = patchFaceRotation[patchFacei];
693 }
694
695 forAll(patchFaceMap, patchFacei)
696 {
697 if (patchFaceMap[patchFacei] != patchFacei)
698 {
699 faceMap[patchFacei + patchStartFace]
700 = patchFaceMap[patchFacei] + patchStartFace;
701 }
702 }
703
704 anyChanged = true;
705 }
706 }
707 }
708
709 Info<< incrIndent << indent << "Faces matched." << endl;
710
711 if (returnReduceOr(anyChanged))
712 {
713 label nReorderedFaces = 0;
714
715 forAll(faceMap, facei)
716 {
717 if (faceMap[facei] != -1)
718 {
719 nReorderedFaces++;
720 }
721 }
722
723 if (nReorderedFaces > 0)
724 {
725 inplaceReorder(faceMap, faces);
726 }
727
728 // Rotate faces (rotation is already in new face indices).
729 label nRotated = 0;
730
731 forAll(rotation, facei)
732 {
733 if (rotation[facei] != 0)
734 {
735 faces[facei] = rotateList(faces[facei], rotation[facei]);
736 nRotated++;
737 }
738 }
739
741 << " faces have been reordered" << nl
742 << indent << returnReduce(nRotated, sumOp<label>())
743 << " faces have been rotated"
745 << decrIndent << decrIndent << endl;
746 }
747}
748
749
751(
752 const word& meshName,
753 const fileName& instance,
755 labelList& boundaryPts,
756 faceList& faces,
757 labelList& owner,
758 labelList& neighbour,
759 const wordList& patchNames,
761 const pointField& cellCentres,
762 bitSet& boundaryFacesToRemove
763) const
764{
765 if (foamyHexMeshControls().objOutput())
766 {
768 (
769 time().path()/word(meshName + ".obj"),
770 points,
771 faces
772 );
773 }
774
775 const label nInternalFaces = patchDicts[0].get<label>("startFace");
776
777 reorderPoints(points, boundaryPts, faces, nInternalFaces);
778
779 if (Pstream::parRun())
780 {
781 reorderProcessorPatches
782 (
783 meshName,
784 instance,
785 points,
786 faces,
789 );
790 }
791
793 Info<< indent << "Constructing mesh" << endl;
794
795 timeCheck("Before fvMesh construction");
796
798 (
800 (
801 meshName,
802 instance,
803 runTime_,
806 ),
807 std::move(points),
808 std::move(faces),
809 std::move(owner),
810 std::move(neighbour)
811 );
812
813 Info<< indent << "Adding patches to mesh" << endl;
814
816
817 label nValidPatches = 0;
818
820 {
821 label totalPatchSize = patchDicts[p].get<label>("nFaces");
822
823 if
824 (
825 patchDicts.set(p)
826 && (
827 patchDicts[p].get<word>("type")
828 == processorPolyPatch::typeName
829 )
830 )
831 {
832 const_cast<dictionary&>(patchDicts[p]).set
833 (
834 "transform",
835 "coincidentFullMatch"
836 );
837
838 // Do not create empty processor patches
839 if (totalPatchSize > 0)
840 {
841 patches[nValidPatches] = new processorPolyPatch
842 (
843 patchNames[p],
844 patchDicts[p],
845 nValidPatches,
847 processorPolyPatch::typeName
848 );
849
850 nValidPatches++;
851 }
852 }
853 else
854 {
855 // Check that the patch is not empty on every processor
856 reduce(totalPatchSize, sumOp<label>());
857
858 if (totalPatchSize > 0)
859 {
860 patches[nValidPatches] = polyPatch::New
861 (
862 patchNames[p],
863 patchDicts[p],
864 nValidPatches,
866 ).ptr();
867
868 nValidPatches++;
869 }
870 }
871 }
872
873 patches.setSize(nValidPatches);
874
876
877
878 // Add zones to the mesh
879 addZones(mesh, cellCentres);
880
881
882 Info<< indent << "Add pointZones" << endl;
883 {
884 label sz = mesh.pointZones().size();
885
886 DynamicList<label> bPts(boundaryPts.size());
887
888 forAll(dualMeshPointTypeNames_, typeI)
889 {
890 const word& znName =
891 dualMeshPointTypeNames_[dualMeshPointType(typeI)];
892
893 forAll(boundaryPts, ptI)
894 {
895 const label& bPtType = boundaryPts[ptI];
896
897 if (bPtType == typeI)
898 {
899 bPts.append(ptI);
900 }
901 }
902
903// syncTools::syncPointList(mesh, bPts, maxEqOp<label>(), -1);
904
906 << "Adding " << bPts.size()
907 << " points of type " << znName
908 << decrIndent << endl;
909
911 (
912 new pointZone
913 (
914 znName,
915 bPts,
916 sz + typeI,
918 )
919 );
920
921 bPts.clear();
922 }
923 }
924
925
926
927 // Add indirectPatchFaces to a face zone
928 Info<< indent << "Adding indirect patch faces set" << endl;
929
931 (
932 mesh,
933 boundaryFacesToRemove,
935 );
936
937 labelList addr(boundaryFacesToRemove.toc());
938
939 faceSet indirectPatchFaces
940 (
941 mesh,
942 "indirectPatchFaces",
943 addr,
945 );
946
947 indirectPatchFaces.sync(mesh);
948
949
951
952 timeCheck("Before fvMesh filtering");
953
954 autoPtr<polyMeshFilter> meshFilter;
955
956 label nInitialBadFaces = 0;
957
958 if (foamyHexMeshControls().filterEdges())
959 {
960 Info<< nl << "Filtering edges on polyMesh" << nl << endl;
961
962 meshFilter.reset(new polyMeshFilter(mesh, boundaryPts));
963
964 // Filter small edges only. This reduces the number of faces so that
965 // the face filtering is sped up.
966 nInitialBadFaces = meshFilter().filterEdges(0);
967 {
968 const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
969
970 polyTopoChange meshMod(newMesh());
971
972 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
973
974 polyMeshFilter::copySets(newMesh(), mesh);
975 }
976 }
977
978 if (foamyHexMeshControls().filterFaces())
979 {
980 labelIOList boundaryPtsIO
981 (
983 (
984 "pointPriority",
985 instance,
986 time(),
989 ),
991 );
992
993 forAll(mesh.points(), ptI)
994 {
995 boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
996 }
997
998
999 Info<< nl << "Filtering faces on polyMesh" << nl << endl;
1000
1001 meshFilter.reset(new polyMeshFilter(mesh, boundaryPtsIO));
1002
1003 meshFilter().filter(nInitialBadFaces);
1004 {
1005 const autoPtr<fvMesh>& newMesh = meshFilter().filteredMesh();
1006
1007 polyTopoChange meshMod(newMesh());
1008
1009 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1010
1011 polyMeshFilter::copySets(newMesh(), mesh);
1012 }
1013 }
1014
1015 timeCheck("After fvMesh filtering");
1016
1017 mesh.setInstance(instance);
1018
1019 if (!mesh.write())
1020 {
1022 << "Failed writing polyMesh."
1023 << exit(FatalError);
1024 }
1025 else
1026 {
1027 Info<< nl << "Written filtered mesh to "
1028 << mesh.polyMesh::instance() << nl
1029 << endl;
1030 }
1031
1032 {
1033 pointScalarField boundaryPtsScalarField
1034 (
1035 IOobject
1036 (
1037 "boundaryPoints_collapsed",
1038 instance,
1039 time(),
1042 ),
1044 dimensionedScalar("min", dimless, scalar(labelMin))
1045 );
1046
1047 labelIOList boundaryPtsIO
1048 (
1049 IOobject
1050 (
1051 "pointPriority",
1052 instance,
1053 time(),
1056 ),
1058 );
1059
1060 forAll(mesh.points(), ptI)
1061 {
1062 boundaryPtsScalarField[ptI] = mesh.pointZones().whichZone(ptI);
1063 boundaryPtsIO[ptI] = mesh.pointZones().whichZone(ptI);
1064 }
1065
1066 boundaryPtsScalarField.write();
1067 boundaryPtsIO.write();
1068 }
1069
1070// writeCellSizes(mesh);
1071
1072// writeCellAlignments(mesh);
1073
1074// writeCellCentres(mesh);
1075
1076 findRemainingProtrusionSet(mesh);
1077}
1078
1079
1081(
1082 const fvMesh& mesh
1083) const
1084{
1085 {
1086 timeCheck("Start writeCellSizes");
1087
1088 Info<< nl << "Create targetCellSize volScalarField" << endl;
1089
1090 volScalarField targetCellSize
1091 (
1092 IOobject
1093 (
1094 "targetCellSize",
1095 mesh.polyMesh::instance(),
1096 mesh,
1099 ),
1100 mesh,
1103 );
1104
1105 scalarField& cellSize = targetCellSize.primitiveFieldRef();
1106
1107 const vectorField& C = mesh.cellCentres();
1108
1109 forAll(cellSize, i)
1110 {
1111 cellSize[i] = cellShapeControls().cellSize(C[i]);
1112 }
1113
1114 // Info<< nl << "Create targetCellVolume volScalarField" << endl;
1115
1116 // volScalarField targetCellVolume
1117 // (
1118 // IOobject
1119 // (
1120 // "targetCellVolume",
1121 // mesh.polyMesh::instance(),
1122 // mesh,
1123 // IOobject::NO_READ,
1124 // IOobject::AUTO_WRITE
1125 // ),
1126 // mesh,
1127 // dimensionedScalar(dimLength, Zero),
1128 // fvPatchFieldBase::zeroGradientType()
1129 // );
1130
1131 // targetCellVolume.primitiveFieldRef() = pow3(cellSize);
1132
1133 // Info<< nl << "Create actualCellVolume volScalarField" << endl;
1134
1135 // volScalarField actualCellVolume
1136 // (
1137 // IOobject
1138 // (
1139 // "actualCellVolume",
1140 // mesh.polyMesh::instance(),
1141 // mesh,
1142 // IOobject::NO_READ,
1143 // IOobject::AUTO_WRITE
1144 // ),
1145 // mesh,
1146 // dimensionedScalar(dimVolume, Zero),
1147 // fvPatchFieldBase::zeroGradientType()
1148 // );
1149
1150 // actualCellVolume.primitiveFieldRef() = mesh.cellVolumes();
1151
1152 // Info<< nl << "Create equivalentCellSize volScalarField" << endl;
1153
1154 // volScalarField equivalentCellSize
1155 // (
1156 // IOobject
1157 // (
1158 // "equivalentCellSize",
1159 // mesh.polyMesh::instance(),
1160 // mesh,
1161 // IOobject::NO_READ,
1162 // IOobject::AUTO_WRITE
1163 // ),
1164 // mesh,
1165 // dimensionedScalar(dimLength, Zero),
1166 // fvPatchFieldBase::zeroGradientType()
1167 // );
1168
1169 // equivalentCellSize.primitiveFieldRef() = pow
1170 // (
1171 // actualCellVolume.primitiveField(),
1172 // 1.0/3.0
1173 // );
1174
1175 targetCellSize.correctBoundaryConditions();
1176 // targetCellVolume.correctBoundaryConditions();
1177 // actualCellVolume.correctBoundaryConditions();
1178 // equivalentCellSize.correctBoundaryConditions();
1179
1180 targetCellSize.write();
1181 // targetCellVolume.write();
1182 // actualCellVolume.write();
1183 // equivalentCellSize.write();
1184 }
1185
1186 // {
1187 // polyMesh tetMesh
1188 // (
1189 // IOobject
1190 // (
1191 // "tetDualMesh",
1192 // runTime_.constant(),
1193 // runTime_,
1194 // IOobject::MUST_READ
1195 // )
1196 // );
1197
1198 // pointMesh ptMesh(tetMesh);
1199
1200 // pointScalarField ptTargetCellSize
1201 // (
1202 // IOobject
1203 // (
1204 // "ptTargetCellSize",
1205 // runTime_.timeName(),
1206 // tetMesh,
1207 // IOobject::NO_READ,
1208 // IOobject::AUTO_WRITE
1209 // ),
1210 // ptMesh,
1211 // dimensionedScalar(dimLength, Zero),
1212 // pointPatchVectorField::calculatedType()
1213 // );
1214
1215 // scalarField& cellSize = ptTargetCellSize.primitiveFieldRef();
1216
1217 // const vectorField& P = tetMesh.points();
1218
1219 // forAll(cellSize, i)
1220 // {
1221 // cellSize[i] = cellShapeControls().cellSize(P[i]);
1222 // }
1223
1224 // ptTargetCellSize.write();
1225 // }
1226}
1227
1228
1230(
1231 const fvMesh& mesh
1232) const
1233{
1234// Info<< nl << "Create cellAlignments volTensorField" << endl;
1235//
1236// volTensorField cellAlignments
1237// (
1238// IOobject
1239// (
1240// "cellAlignments",
1241// mesh.polyMesh::instance(),
1242// mesh,
1243// IOobject::NO_READ,
1244// IOobject::NO_WRITE,
1245// IOobject::NO_REGISTER
1246// ),
1247// mesh,
1248// tensor::I,
1249// fvPatchFieldBase::zeroGradientType()
1250// );
1251//
1252// tensorField& cellAlignment = cellAlignments.primitiveFieldRef();
1253//
1254// const vectorField& C = mesh.cellCentres();
1255//
1256// vectorField xDir(cellAlignment.size());
1257// vectorField yDir(cellAlignment.size());
1258// vectorField zDir(cellAlignment.size());
1259//
1260// forAll(cellAlignment, i)
1261// {
1262// cellAlignment[i] = cellShapeControls().cellAlignment(C[i]);
1263// xDir[i] = cellAlignment[i] & vector(1, 0, 0);
1264// yDir[i] = cellAlignment[i] & vector(0, 1, 0);
1265// zDir[i] = cellAlignment[i] & vector(0, 0, 1);
1266// }
1267//
1268// OFstream xStr("xDir.obj");
1269// OFstream yStr("yDir.obj");
1270// OFstream zStr("zDir.obj");
1271//
1272// forAll(xDir, i)
1273// {
1274// meshTools::writeOBJ(xStr, C[i], C[i] + xDir[i]);
1275// meshTools::writeOBJ(yStr, C[i], C[i] + yDir[i]);
1276// meshTools::writeOBJ(zStr, C[i], C[i] + zDir[i]);
1277// }
1278//
1279// cellAlignments.correctBoundaryConditions();
1280//
1281// cellAlignments.write();
1282}
1283
1284
1286(
1287 const fvMesh& mesh
1288) const
1289{
1290 Info<< "Writing components of cellCentre positions to volScalarFields"
1291 << " ccx, ccy, ccz in " << runTime_.timeName() << endl;
1292
1293 for (direction i=0; i<vector::nComponents; i++)
1294 {
1295 volScalarField cci
1296 (
1297 IOobject
1298 (
1299 "cc" + word(vector::componentNames[i]),
1300 runTime_.timeName(),
1301 mesh,
1304 ),
1305 mesh.C().component(i)
1306 );
1307
1308 cci.write();
1309 }
1310}
1311
1312
1314(
1315 const polyMesh& mesh
1316) const
1317{
1318 timeCheck("Start findRemainingProtrusionSet");
1319
1321
1322 labelHashSet protrudingBoundaryPoints;
1323
1324 forAll(patches, patchi)
1325 {
1326 const polyPatch& patch = patches[patchi];
1327
1328 forAll(patch.localPoints(), pLPI)
1329 {
1330 label meshPtI = patch.meshPoints()[pLPI];
1331
1332 const Foam::point& pt = patch.localPoints()[pLPI];
1333
1334 if
1335 (
1336 geometryToConformTo_.wellOutside
1337 (
1338 pt,
1339 sqr(targetCellSize(pt))
1340 )
1341 )
1342 {
1343 protrudingBoundaryPoints.insert(meshPtI);
1344 }
1345 }
1346 }
1347
1348 cellSet protrudingCells
1349 (
1350 mesh,
1351 "foamyHexMesh_remainingProtrusions",
1352 mesh.nCells()/1000
1353 );
1354
1355 for (const label pointi : protrudingBoundaryPoints)
1356 {
1357 const labelList& pCells = mesh.pointCells()[pointi];
1358 protrudingCells.insert(pCells);
1359 }
1360
1361 const label protrudingCellsSize =
1362 returnReduce(protrudingCells.size(), sumOp<label>());
1363
1364 if (foamyHexMeshControls().objOutput() && protrudingCellsSize)
1365 {
1366 Info<< nl << "Found " << protrudingCellsSize
1367 << " cells protruding from the surface, writing cellSet "
1368 << protrudingCells.name()
1369 << endl;
1370
1371 protrudingCells.write();
1372 }
1373
1374 return protrudingCells;
1375}
1376
1377
1378void Foam::conformalVoronoiMesh::writePointPairs
1379(
1380 const fileName& fName
1381) const
1382{
1383 OBJstream os(fName);
1384
1385 for
1386 (
1387 Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1388 eit != finite_edges_end();
1389 ++eit
1390 )
1391 {
1392 Cell_handle c = eit->first;
1393 Vertex_handle vA = c->vertex(eit->second);
1394 Vertex_handle vB = c->vertex(eit->third);
1395
1396 if (ptPairs_.isPointPair(vA, vB))
1397 {
1398 os.writeLine(topoint(vA->point()), topoint(vB->point()));
1399 }
1400 }
1401}
1402
1403
1404// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Various functions to operate on Lists.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Graphite solid properties.
Definition C.H:49
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ 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
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
Print a table in parallel, e.g.;.
Definition PrintTable.H:67
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition PtrList.H:364
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
const T * get(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrListI.H:134
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const char *const componentNames[]
static constexpr direction nComponents
Number of components in this vector space.
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition ZoneMesh.C:410
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 bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A collection of cell labels.
Definition cellSet.H:50
void writeMesh(const fileName &instance)
Prepare data and call writeMesh for polyMesh and.
labelHashSet findRemainingProtrusionSet(const polyMesh &mesh) const
Find the cellSet of the boundary cells which have points that.
const Time & time() const
Return the Time object.
void writeCellSizes(const fvMesh &mesh) const
Calculate and write a field of the target cell size,.
void writeCellCentres(const fvMesh &mesh) const
Calculate and write the cell centres.
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
void writeCellAlignments(const fvMesh &mesh) const
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
double cpuTimeIncrement() const
Return CPU time [seconds] since last call to cpuTimeIncrement(), resetCpuTimeIncrement().
double elapsedCpuTime() const
Return CPU time [seconds] from the start.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A list of face labels.
Definition faceSet.H:50
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const volVectorField & C() const
Return cell centres as volVectorField.
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition fvMesh.C:628
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition fvMesh.C:1068
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Memory usage information for the current process, and the system memory that is free.
Definition memInfo.H:59
A subset of mesh points.
Definition pointZone.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Remove the edges and faces of a polyMesh whilst satisfying the given mesh quality criteria.
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
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
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.
Direct mesh changes based on v1.3 polyTopoChange syntax.
const labelListList & pointCells() const
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
Neighbour processor patch.
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
A class for handling words, derived from Foam::string.
Definition word.H:66
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
const pointField & points
word timeName
Definition getTimeIndex.H:3
void writeInternalDelaunayVertices(const fileName &instance, const Triangulation &t)
Write the internal Delaunay vertices of the tessellation as a.
void writeObjMesh(const fileName &fName, const pointField &points, const faceList &faces)
Write an OBJ mesh consisting of points and faces.
const dimensionedScalar c
Speed of light in a vacuum.
const std::string patch
OpenFOAM patch number as a std::string.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
constexpr label labelMin
Definition label.H:54
GeometricField< scalar, fvPatchField, volMesh > volScalarField
pointFromPoint topoint(const Point &P)
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
ListType rotateList(const ListType &list, const label n)
Rotate a list by n places.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
static void check(const int retVal, const char *what)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
HashTable< label, labelPair, Foam::Hash< labelPair > > labelPairLookup
This is a Map of a labelPair to a label. Used for e.g. for face1, face2 to shared edge....
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
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).
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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...
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
PtrList< dictionary > patchDicts
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299