Loading...
Searching...
No Matches
meshDualiser.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) 2020-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "meshDualiser.H"
30#include "meshTools.H"
31#include "polyMesh.H"
32#include "polyTopoChange.H"
33#include "mapPolyMesh.H"
34#include "edgeFaceCirculator.H"
35#include "mergePoints.H"
36#include "DynamicList.H"
37#include "OFstream.H"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::meshDualiser::checkPolyTopoChange(const polyTopoChange& meshMod)
50{
51 // Assume no removed points
52 pointField points(meshMod.points().size());
53 forAll(meshMod.points(), i)
54 {
55 points[i] = meshMod.points()[i];
56 }
57
58 labelList oldToNew;
59 label nUnique = Foam::mergePoints
60 (
61 points,
62 1e-6,
63 false,
64 oldToNew
65 );
66
67 if (nUnique < points.size())
68 {
69 CompactListList<label> newToOld
70 (
71 invertOneToManyCompact(nUnique, oldToNew)
72 );
73
74 forAll(newToOld, newi)
75 {
76 if (newToOld[newi].size() != 1)
77 {
79 << "duplicate verts:" << newToOld[newi]
80 << " coords:"
81 << UIndirectList<point>(points, newToOld[newi])
82 << abort(FatalError);
83 }
84 }
85 }
86}
87
88
89// Dump state so far.
90void Foam::meshDualiser::dumpPolyTopoChange
91(
92 const polyTopoChange& meshMod,
93 const fileName& prefix
94)
95{
96 OFstream str1(prefix + "Faces.obj");
97 OFstream str2(prefix + "Edges.obj");
98
99 Info<< "Dumping current polyTopoChange. Faces to " << str1.name()
100 << " , points and edges to " << str2.name() << endl;
101
102 for (const auto& p : meshMod.points())
103 {
104 meshTools::writeOBJ(str1, p);
105 meshTools::writeOBJ(str2, p);
106 }
107
108 for (const face& f : meshMod.faces())
109 {
110 str1<< 'f';
111 forAll(f, fp)
112 {
113 str1<< ' ' << f[fp]+1;
114 }
115 str1<< nl;
116
117 str2<< 'l';
118 forAll(f, fp)
119 {
120 str2<< ' ' << f[fp]+1;
121 }
122 str2<< ' ' << f[0]+1 << nl;
123 }
124}
125
126
127Foam::label Foam::meshDualiser::findDualCell
128(
129 const label celli,
130 const label pointi
131) const
132{
133 const labelList& dualCells = pointToDualCells_[pointi];
134
135 if (dualCells.size() == 1)
136 {
137 return dualCells[0];
138 }
139 else
140 {
141 label index = mesh_.pointCells()[pointi].find(celli);
142
143 return dualCells[index];
144 }
145}
146
147
148void Foam::meshDualiser::generateDualBoundaryEdges
149(
150 const bitSet& isBoundaryEdge,
151 const label pointi,
152 polyTopoChange& meshMod
153)
154{
155 const labelList& pEdges = mesh_.pointEdges()[pointi];
156
157 forAll(pEdges, pEdgeI)
158 {
159 label edgeI = pEdges[pEdgeI];
160
161 if (edgeToDualPoint_[edgeI] == -1 && isBoundaryEdge.test(edgeI))
162 {
163 const edge& e = mesh_.edges()[edgeI];
164
165 edgeToDualPoint_[edgeI] = meshMod.addPoint
166 (
167 e.centre(mesh_.points()),
168 pointi, // masterPoint
169 -1, // zoneID
170 true // inCell
171 );
172 }
173 }
174}
175
176
177// Return true if point on face has same dual cells on both owner and neighbour
178// sides.
179bool Foam::meshDualiser::sameDualCell
180(
181 const label facei,
182 const label pointi
183) const
184{
185 if (!mesh_.isInternalFace(facei))
186 {
188 << "face:" << facei << " is not internal face."
189 << abort(FatalError);
190 }
191
192 label own = mesh_.faceOwner()[facei];
193 label nei = mesh_.faceNeighbour()[facei];
194
195 return findDualCell(own, pointi) == findDualCell(nei, pointi);
196}
197
198
199Foam::label Foam::meshDualiser::addInternalFace
200(
201 const label masterPointi,
202 const label masterEdgeI,
203 const label masterFacei,
204
205 const bool edgeOrder,
206 const label dualCell0,
207 const label dualCell1,
208 const labelUList& verts,
209 polyTopoChange& meshMod
210) const
211{
212 face newFace(verts);
213
214 if (edgeOrder != (dualCell0 < dualCell1))
215 {
216 reverse(newFace);
217 }
218
219 if (debug)
220 {
221 pointField facePoints(meshMod.points(), newFace);
222
223 labelList oldToNew;
224 label nUnique = Foam::mergePoints
225 (
226 facePoints,
227 1e-6,
228 false,
229 oldToNew
230 );
231
232 if (nUnique < facePoints.size())
233 {
235 << "verts:" << verts << " newFace:" << newFace
236 << " face points:" << facePoints
237 << abort(FatalError);
238 }
239 }
240
241
242 label zoneID = -1;
243 bool zoneFlip = false;
244 if (masterFacei != -1)
245 {
246 zoneID = mesh_.faceZones().whichZone(masterFacei);
247
248 if (zoneID != -1)
249 {
250 const faceZone& fZone = mesh_.faceZones()[zoneID];
251
252 zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
253 }
254 }
255
256 label dualFacei;
257
258 if (dualCell0 < dualCell1)
259 {
260 dualFacei = meshMod.addFace
261 (
262 newFace,
263 dualCell0, // own
264 dualCell1, // nei
265 masterPointi, // masterPointID
266 masterEdgeI, // masterEdgeID
267 masterFacei, // masterFaceID
268 false, // flipFaceFlux
269 -1, // patchID
270 zoneID, // zoneID
271 zoneFlip // zoneFlip
272 );
273
274 //pointField dualPoints(meshMod.points());
275 //const vector n(newFace.unitNormal(dualPoints));
276 //
277 //Pout<< "Generated internal dualFace:" << dualFacei
278 // << " verts:" << newFace
279 // << " points:" << UIndirectList<point>(meshMod.points(), newFace)
280 // << " n:" << n
281 // << " between dualowner:" << dualCell0
282 // << " dualneighbour:" << dualCell1
283 // << endl;
284 }
285 else
286 {
287 dualFacei = meshMod.addFace
288 (
289 newFace,
290 dualCell1, // own
291 dualCell0, // nei
292 masterPointi, // masterPointID
293 masterEdgeI, // masterEdgeID
294 masterFacei, // masterFaceID
295 false, // flipFaceFlux
296 -1, // patchID
297 zoneID, // zoneID
298 zoneFlip // zoneFlip
299 );
300
301 //pointField dualPoints(meshMod.points());
302 //const vector n(newFace.unitNormal(dualPoints));
303 //
304 //Pout<< "Generated internal dualFace:" << dualFacei
305 // << " verts:" << newFace
306 // << " points:" << UIndirectList<point>(meshMod.points(), newFace)
307 // << " n:" << n
308 // << " between dualowner:" << dualCell1
309 // << " dualneighbour:" << dualCell0
310 // << endl;
311 }
312 return dualFacei;
313}
314
315
316Foam::label Foam::meshDualiser::addBoundaryFace
317(
318 const label masterPointi,
319 const label masterEdgeI,
320 const label masterFacei,
321
322 const label dualCelli,
323 const label patchi,
324 const labelUList& verts,
325 polyTopoChange& meshMod
326) const
327{
328 face newFace(verts);
329
330 label zoneID = -1;
331 bool zoneFlip = false;
332 if (masterFacei != -1)
333 {
334 zoneID = mesh_.faceZones().whichZone(masterFacei);
335
336 if (zoneID != -1)
337 {
338 const faceZone& fZone = mesh_.faceZones()[zoneID];
339
340 zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
341 }
342 }
343
344 label dualFacei = meshMod.addFace
345 (
346 newFace,
347 dualCelli, // own
348 -1, // nei
349 masterPointi, // masterPointID
350 masterEdgeI, // masterEdgeID
351 masterFacei, // masterFaceID
352 false, // flipFaceFlux
353 patchi, // patchID
354 zoneID, // zoneID
355 zoneFlip // zoneFlip
356 );
357
358 //pointField dualPoints(meshMod.points());
359 //const vector n(newFace.unitNormal(dualPoints));
360 //
361 //Pout<< "Generated boundary dualFace:" << dualFacei
362 // << " verts:" << newFace
363 // << " points:" << UIndirectList<point>(meshMod.points(), newFace)
364 // << " n:" << n
365 // << " on dualowner:" << dualCelli
366 // << endl;
367 return dualFacei;
368}
369
370
371// Walks around edgeI.
372// splitFace=true : creates multiple faces
373// splitFace=false: creates single face if same dual cells on both sides,
374// multiple faces otherwise.
375void Foam::meshDualiser::createFacesAroundEdge
376(
377 const bool splitFace,
378 const bitSet& isBoundaryEdge,
379 const label edgeI,
380 const label startFacei,
381 polyTopoChange& meshMod,
382 boolList& doneEFaces
383) const
384{
385 const edge& e = mesh_.edges()[edgeI];
386 const labelList& eFaces = mesh_.edgeFaces()[edgeI];
387
389 (
390 mesh_.faces()[startFacei],
391 e[0],
392 e[1]
393 );
394
396 (
397 mesh_,
398 startFacei, // face
399 true, // ownerSide
400 fp, // fp
401 isBoundaryEdge.test(edgeI) // isBoundaryEdge
402 );
403 ie.setCanonical();
404
405 bool edgeOrder = ie.sameOrder(e[0], e[1]);
406 label startFaceLabel = ie.faceLabel();
407
408 //Pout<< "At edge:" << edgeI << " verts:" << e
409 // << " points:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
410 // << " started walking at face:" << ie.faceLabel()
411 // << " verts:" << mesh_.faces()[ie.faceLabel()]
412 // << " edgeOrder:" << edgeOrder
413 // << " in direction of cell:" << ie.cellLabel()
414 // << endl;
415
416 // Walk and collect face.
417 DynamicList<label> verts(100);
418
419 if (edgeToDualPoint_[edgeI] != -1)
420 {
421 verts.append(edgeToDualPoint_[edgeI]);
422 }
423 if (faceToDualPoint_[ie.faceLabel()] != -1)
424 {
425 doneEFaces[eFaces.find(ie.faceLabel())] = true;
426 verts.append(faceToDualPoint_[ie.faceLabel()]);
427 }
428 if (cellToDualPoint_[ie.cellLabel()] != -1)
429 {
430 verts.append(cellToDualPoint_[ie.cellLabel()]);
431 }
432
433 label currentDualCell0 = findDualCell(ie.cellLabel(), e[0]);
434 label currentDualCell1 = findDualCell(ie.cellLabel(), e[1]);
435
436 ++ie;
437
438 while (true)
439 {
440 label facei = ie.faceLabel();
441
442 // Mark face as visited.
443 doneEFaces[eFaces.find(facei)] = true;
444
445 if (faceToDualPoint_[facei] != -1)
446 {
447 verts.append(faceToDualPoint_[facei]);
448 }
449
450 label celli = ie.cellLabel();
451
452 if (celli == -1)
453 {
454 // At ending boundary face. We've stored the face point above
455 // so this is the whole face.
456 break;
457 }
458
459
460 label dualCell0 = findDualCell(celli, e[0]);
461 label dualCell1 = findDualCell(celli, e[1]);
462
463 // Generate face. (always if splitFace=true; only if needed to
464 // separate cells otherwise)
465 if
466 (
467 splitFace
468 || (
469 dualCell0 != currentDualCell0
470 || dualCell1 != currentDualCell1
471 )
472 )
473 {
474 // Close current face.
475 addInternalFace
476 (
477 -1, // masterPointi
478 edgeI, // masterEdgeI
479 -1, // masterFacei
480 edgeOrder,
481 currentDualCell0,
482 currentDualCell1,
483 verts.shrink(),
484 meshMod
485 );
486
487 // Restart
488 currentDualCell0 = dualCell0;
489 currentDualCell1 = dualCell1;
490
491 verts.clear();
492 if (edgeToDualPoint_[edgeI] != -1)
493 {
494 verts.append(edgeToDualPoint_[edgeI]);
495 }
496 if (faceToDualPoint_[facei] != -1)
497 {
498 verts.append(faceToDualPoint_[facei]);
499 }
500 }
501
502 if (cellToDualPoint_[celli] != -1)
503 {
504 verts.append(cellToDualPoint_[celli]);
505 }
506
507 ++ie;
508
509 if (ie == ie.end())
510 {
511 // Back at start face (for internal edge only). See if this needs
512 // adding.
513 if (!isBoundaryEdge.test(edgeI))
514 {
515 label startDual = faceToDualPoint_[startFaceLabel];
516
517 if (startDual != -1)
518 {
519 verts.push_uniq(startDual);
520 }
521 }
522 break;
523 }
524 }
525
526 verts.shrink();
527 addInternalFace
528 (
529 -1, // masterPointi
530 edgeI, // masterEdgeI
531 -1, // masterFacei
532 edgeOrder,
533 currentDualCell0,
534 currentDualCell1,
535 verts,
536 meshMod
537 );
538}
539
540
541// Walks around circumference of facei. Creates single face. Gets given
542// starting (feature) edge to start from. Returns ending edge. (all edges
543// in form of index in faceEdges)
544void Foam::meshDualiser::createFaceFromInternalFace
545(
546 const label facei,
547 label& fp,
548 polyTopoChange& meshMod
549) const
550{
551 const face& f = mesh_.faces()[facei];
552 const labelList& fEdges = mesh_.faceEdges()[facei];
553 label own = mesh_.faceOwner()[facei];
554 label nei = mesh_.faceNeighbour()[facei];
555
556 //Pout<< "createFaceFromInternalFace : At face:" << facei
557 // << " verts:" << f
558 // << " points:" << UIndirectList<point>(mesh_.points(), f)
559 // << " started walking at edge:" << fEdges[fp]
560 // << " verts:" << mesh_.edges()[fEdges[fp]]
561 // << endl;
562
563
564 // Walk and collect face.
565 DynamicList<label> verts(100);
566
567 verts.append(faceToDualPoint_[facei]);
568 verts.append(edgeToDualPoint_[fEdges[fp]]);
569
570 // Step to vertex after edge mid
571 fp = f.fcIndex(fp);
572
573 // Get cells on either side of face at that point
574 label currentDualCell0 = findDualCell(own, f[fp]);
575 label currentDualCell1 = findDualCell(nei, f[fp]);
576
577 forAll(f, i)
578 {
579 // Check vertex
580 if (pointToDualPoint_[f[fp]] != -1)
581 {
582 verts.append(pointToDualPoint_[f[fp]]);
583 }
584
585 // Edge between fp and fp+1
586 label edgeI = fEdges[fp];
587
588 if (edgeToDualPoint_[edgeI] != -1)
589 {
590 verts.append(edgeToDualPoint_[edgeI]);
591 }
592
593 // Next vertex on edge
594 label nextFp = f.fcIndex(fp);
595
596 // Get dual cells on nextFp to check whether face needs closing.
597 label dualCell0 = findDualCell(own, f[nextFp]);
598 label dualCell1 = findDualCell(nei, f[nextFp]);
599
600 if (dualCell0 != currentDualCell0 || dualCell1 != currentDualCell1)
601 {
602 // Check: make sure that there is a midpoint on the edge.
603 if (edgeToDualPoint_[edgeI] == -1)
604 {
606 << "face:" << facei << " verts:" << f
607 << " points:" << UIndirectList<point>(mesh_.points(), f)
608 << " no feature edge between " << f[fp]
609 << " and " << f[nextFp] << " although have different"
610 << " dual cells." << endl
611 << "point " << f[fp] << " has dual cells "
612 << currentDualCell0 << " and " << currentDualCell1
613 << " ; point "<< f[nextFp] << " has dual cells "
614 << dualCell0 << " and " << dualCell1
615 << abort(FatalError);
616 }
617
618
619 // Close current face.
620 verts.shrink();
621 addInternalFace
622 (
623 -1, // masterPointi
624 -1, // masterEdgeI
625 facei, // masterFacei
626 true, // edgeOrder,
627 currentDualCell0,
628 currentDualCell1,
629 verts,
630 meshMod
631 );
632 break;
633 }
634
635 fp = nextFp;
636 }
637}
638
639
640// Given a point on a face converts the faces around the point.
641// (pointFaces()). Gets starting face and marks off visited faces in donePFaces.
642void Foam::meshDualiser::createFacesAroundBoundaryPoint
643(
644 const label patchi,
645 const label patchPointi,
646 const label startFacei,
647 polyTopoChange& meshMod,
648 boolList& donePFaces // pFaces visited
649) const
650{
651 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
652 const polyPatch& pp = patches[patchi];
653 const labelList& pFaces = pp.pointFaces()[patchPointi];
654 const labelList& own = mesh_.faceOwner();
655
656 label pointi = pp.meshPoints()[patchPointi];
657
658 if (pointToDualPoint_[pointi] == -1)
659 {
660 // Not a feature point. Loop over all connected
661 // pointFaces.
662
663 // Starting face
664 label facei = startFacei;
665
666 DynamicList<label> verts(4);
667
668 while (true)
669 {
670 label index = pFaces.find(facei-pp.start());
671
672 // Has face been visited already?
673 if (donePFaces[index])
674 {
675 break;
676 }
677 donePFaces[index] = true;
678
679 // Insert face centre
680 verts.append(faceToDualPoint_[facei]);
681
682 label dualCelli = findDualCell(own[facei], pointi);
683
684 // Get the edge before the patchPointi
685 const face& f = mesh_.faces()[facei];
686 label fp = f.find(pointi);
687 label prevFp = f.rcIndex(fp);
688 label edgeI = mesh_.faceEdges()[facei][prevFp];
689
690 if (edgeToDualPoint_[edgeI] != -1)
691 {
692 verts.append(edgeToDualPoint_[edgeI]);
693 }
694
695 // Get next boundary face (whilst staying on edge).
697 (
698 mesh_,
699 facei,
700 true, // ownerSide
701 prevFp, // index of edge in face
702 true // isBoundaryEdge
703 );
704
705 do
706 {
707 ++circ;
708 }
709 while (mesh_.isInternalFace(circ.faceLabel()));
710
711 // Step to next face
712 facei = circ.faceLabel();
713
714 if (facei < pp.start() || facei >= pp.start()+pp.size())
715 {
717 << "Walked from face on patch:" << patchi
718 << " to face:" << facei
719 << " fc:" << mesh_.faceCentres()[facei]
720 << " on patch:" << patches.whichPatch(facei)
721 << abort(FatalError);
722 }
723
724 // Check if different cell.
725 if (dualCelli != findDualCell(own[facei], pointi))
726 {
728 << "Different dual cells but no feature edge"
729 << " inbetween point:" << pointi
730 << " coord:" << mesh_.points()[pointi]
731 << abort(FatalError);
732 }
733 }
734
735 verts.shrink();
736
737 label dualCelli = findDualCell(own[facei], pointi);
738
739 //Bit dodgy: create dualface from the last face (instead of from
740 // the central point). This will also use the original faceZone to
741 // put the new face (which might span multiple original faces) in.
742
743 addBoundaryFace
744 (
745 //pointi, // masterPointi
746 -1, // masterPointi
747 -1, // masterEdgeI
748 facei, // masterFacei
749 dualCelli,
750 patchi,
751 verts,
752 meshMod
753 );
754 }
755 else
756 {
757 label facei = startFacei;
758
759 // Storage for face
760 DynamicList<label> verts(mesh_.faces()[facei].size());
761
762 // Starting point.
763 verts.append(pointToDualPoint_[pointi]);
764
765 // Find edge between pointi and next point on face.
766 const labelList& fEdges = mesh_.faceEdges()[facei];
767 label nextEdgeI = fEdges[mesh_.faces()[facei].find(pointi)];
768 if (edgeToDualPoint_[nextEdgeI] != -1)
769 {
770 verts.append(edgeToDualPoint_[nextEdgeI]);
771 }
772
773 do
774 {
775 label index = pFaces.find(facei-pp.start());
776
777 // Has face been visited already?
778 if (donePFaces[index])
779 {
780 break;
781 }
782 donePFaces[index] = true;
783
784 // Face centre
785 verts.append(faceToDualPoint_[facei]);
786
787 // Find edge before pointi on facei
788 const labelList& fEdges = mesh_.faceEdges()[facei];
789 const face& f = mesh_.faces()[facei];
790 label prevFp = f.rcIndex(f.find(pointi));
791 label edgeI = fEdges[prevFp];
792
793 if (edgeToDualPoint_[edgeI] != -1)
794 {
795 // Feature edge. Close any face so far. Note: uses face to
796 // create dualFace from. Could use pointi instead.
797 verts.append(edgeToDualPoint_[edgeI]);
798 addBoundaryFace
799 (
800 -1, // masterPointi
801 -1, // masterEdgeI
802 facei, // masterFacei
803 findDualCell(own[facei], pointi),
804 patchi,
805 verts.shrink(),
806 meshMod
807 );
808 verts.clear();
809
810 verts.append(pointToDualPoint_[pointi]);
811 verts.append(edgeToDualPoint_[edgeI]);
812 }
813
814 // Cross edgeI to next boundary face
816 (
817 mesh_,
818 facei,
819 true, // ownerSide
820 prevFp, // index of edge in face
821 true // isBoundaryEdge
822 );
823
824 do
825 {
826 ++circ;
827 }
828 while (mesh_.isInternalFace(circ.faceLabel()));
829
830 // Step to next face. Quit if not on same patch.
831 facei = circ.faceLabel();
832 }
833 while
834 (
835 facei != startFacei
836 && facei >= pp.start()
837 && facei < pp.start()+pp.size()
838 );
839
840 if (verts.size() > 2)
841 {
842 // Note: face created from face, not from pointi
843 addBoundaryFace
844 (
845 -1, // masterPointi
846 -1, // masterEdgeI
847 startFacei, // masterFacei
848 findDualCell(own[facei], pointi),
849 patchi,
850 verts.shrink(),
851 meshMod
852 );
853 }
854 }
855}
856
857
858// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
859
860Foam::meshDualiser::meshDualiser(const polyMesh& mesh)
861:
862 mesh_(mesh),
863 pointToDualCells_(mesh_.nPoints()),
864 pointToDualPoint_(mesh_.nPoints(), -1),
865 cellToDualPoint_(mesh_.nCells()),
866 faceToDualPoint_(mesh_.nFaces(), -1),
867 edgeToDualPoint_(mesh_.nEdges(), -1)
868{}
869
870
871// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
872
874(
875 const bool splitFace,
876 const labelList& featureFaces,
877 const labelList& featureEdges,
878 const labelList& singleCellFeaturePoints,
879 const labelList& multiCellFeaturePoints,
880 polyTopoChange& meshMod
881)
882{
883 const labelList& own = mesh_.faceOwner();
884 const labelList& nei = mesh_.faceNeighbour();
885 const vectorField& cellCentres = mesh_.cellCentres();
886
887 // Mark boundary edges and points.
888 // (Note: in 1.4.2 we can use the built-in mesh point ordering
889 // facility instead)
890 bitSet isBoundaryEdge(mesh_.nEdges());
891 for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
892 {
893 const labelList& fEdges = mesh_.faceEdges()[facei];
894
895 isBoundaryEdge.set(fEdges);
896 }
897
898
899 if (splitFace)
900 {
901 // This is a special mode where whenever we are walking around an edge
902 // every area through a cell becomes a separate dualface. So two
903 // dual cells will probably have more than one dualface between them!
904 // This mode implies that
905 // - all faces have to be feature faces since there has to be a
906 // dualpoint at the face centre.
907 // - all edges have to be feature edges ,,
908 boolList featureFaceSet(mesh_.nFaces(), false);
909 forAll(featureFaces, i)
910 {
911 featureFaceSet[featureFaces[i]] = true;
912 }
913 label facei = featureFaceSet.find(false);
914
915 if (facei != -1)
916 {
918 << "In split-face-mode (splitFace=true) but not all faces"
919 << " marked as feature faces." << endl
920 << "First conflicting face:" << facei
921 << " centre:" << mesh_.faceCentres()[facei]
922 << abort(FatalError);
923 }
924
925 boolList featureEdgeSet(mesh_.nEdges(), false);
926 forAll(featureEdges, i)
927 {
928 featureEdgeSet[featureEdges[i]] = true;
929 }
930 label edgeI = featureEdgeSet.find(false);
931
932 if (edgeI != -1)
933 {
934 const edge& e = mesh_.edges()[edgeI];
936 << "In split-face-mode (splitFace=true) but not all edges"
937 << " marked as feature edges." << endl
938 << "First conflicting edge:" << edgeI
939 << " verts:" << e
940 << " coords:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
941 << abort(FatalError);
942 }
943 }
944 else
945 {
946 // Check that all boundary faces are feature faces.
947
948 boolList featureFaceSet(mesh_.nFaces(), false);
949 forAll(featureFaces, i)
950 {
951 featureFaceSet[featureFaces[i]] = true;
952 }
953 for
954 (
955 label facei = mesh_.nInternalFaces();
956 facei < mesh_.nFaces();
957 facei++
958 )
959 {
960 if (!featureFaceSet[facei])
961 {
963 << "Not all boundary faces marked as feature faces."
964 << endl
965 << "First conflicting face:" << facei
966 << " centre:" << mesh_.faceCentres()[facei]
967 << abort(FatalError);
968 }
969 }
970 }
971
972
973
974
975 // Start creating cells, points, and faces (in that order)
976
977
978 // 1. Mark which cells to create
979 // Mostly every point becomes one cell but sometimes (for feature points)
980 // all cells surrounding a feature point become cells. Also a non-manifold
981 // point can create two cells! So a dual cell is uniquely defined by a
982 // mesh point + cell (as in pointCells index)
983
984 // 2. Mark which face centres to create
985
986 // 3. Internal faces can now consist of
987 // - only cell centres of walk around edge
988 // - cell centres + face centres of walk around edge
989 // - same but now other side is not a single cell
990
991 // 4. Boundary faces (or internal faces between cell zones!) now consist of
992 // - walk around boundary point.
993
994
995
996 autoPtr<OFstream> dualCcStr;
997 if (debug)
998 {
999 dualCcStr.reset(new OFstream("dualCc.obj"));
1000 Pout<< "Dumping centres of dual cells to " << dualCcStr().name()
1001 << endl;
1002 }
1003
1004
1005 // Dual cells (from points)
1006 // ~~~~~~~~~~~~~~~~~~~~~~~~
1007
1008 // pointToDualCells_[pointi]
1009 // - single entry : all cells surrounding point all become the same
1010 // cell.
1011 // - multiple entries: in order of pointCells.
1012
1013
1014 // feature points that become single cell
1015 forAll(singleCellFeaturePoints, i)
1016 {
1017 label pointi = singleCellFeaturePoints[i];
1018
1019 pointToDualPoint_[pointi] = meshMod.addPoint
1020 (
1021 mesh_.points()[pointi],
1022 pointi, // masterPoint
1023 mesh_.pointZones().whichZone(pointi), // zoneID
1024 true // inCell
1025 );
1026
1027 // Generate single cell
1028 pointToDualCells_[pointi].setSize(1);
1029 pointToDualCells_[pointi][0] = meshMod.addCell
1030 (
1031 pointi, //masterPointID,
1032 -1, //masterEdgeID,
1033 -1, //masterFaceID,
1034 -1, //masterCellID,
1035 -1 //zoneID
1036 );
1037 if (dualCcStr)
1038 {
1039 meshTools::writeOBJ(*dualCcStr, mesh_.points()[pointi]);
1040 }
1041 }
1042
1043 // feature points that become multiple cells
1044 forAll(multiCellFeaturePoints, i)
1045 {
1046 label pointi = multiCellFeaturePoints[i];
1047
1048 if (pointToDualCells_[pointi].size() > 0)
1049 {
1051 << "Point " << pointi << " at:" << mesh_.points()[pointi]
1052 << " is both in singleCellFeaturePoints"
1053 << " and multiCellFeaturePoints."
1054 << abort(FatalError);
1055 }
1056
1057 pointToDualPoint_[pointi] = meshMod.addPoint
1058 (
1059 mesh_.points()[pointi],
1060 pointi, // masterPoint
1061 mesh_.pointZones().whichZone(pointi), // zoneID
1062 true // inCell
1063 );
1064
1065 // Create dualcell for every cell connected to dual point
1066
1067 const labelList& pCells = mesh_.pointCells()[pointi];
1068
1069 pointToDualCells_[pointi].setSize(pCells.size());
1070
1071 forAll(pCells, pCelli)
1072 {
1073 pointToDualCells_[pointi][pCelli] = meshMod.addCell
1074 (
1075 pointi, //masterPointID
1076 -1, //masterEdgeID
1077 -1, //masterFaceID
1078 -1, //masterCellID
1079 mesh_.cellZones().whichZone(pCells[pCelli]) //zoneID
1080 );
1081 if (dualCcStr)
1082 {
1084 (
1085 *dualCcStr,
1086 0.5*(mesh_.points()[pointi]+cellCentres[pCells[pCelli]])
1087 );
1088 }
1089 }
1090 }
1091 // Normal points
1092 forAll(mesh_.points(), pointi)
1093 {
1094 if (pointToDualCells_[pointi].empty())
1095 {
1096 pointToDualCells_[pointi].setSize(1);
1097 pointToDualCells_[pointi][0] = meshMod.addCell
1098 (
1099 pointi, //masterPointID,
1100 -1, //masterEdgeID,
1101 -1, //masterFaceID,
1102 -1, //masterCellID,
1103 -1 //zoneID
1104 );
1105
1106 if (dualCcStr)
1107 {
1108 meshTools::writeOBJ(*dualCcStr, mesh_.points()[pointi]);
1109 }
1110 }
1111 }
1112
1113
1114 // Dual points (from cell centres, feature faces, feature edges)
1115 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1116
1117 forAll(cellToDualPoint_, celli)
1118 {
1119 cellToDualPoint_[celli] = meshMod.addPoint
1120 (
1121 cellCentres[celli],
1122 mesh_.faces()[mesh_.cells()[celli][0]][0], // masterPoint
1123 -1, // zoneID
1124 true // inCell
1125 );
1126 }
1127
1128 // From face to dual point
1129
1130 forAll(featureFaces, i)
1131 {
1132 label facei = featureFaces[i];
1133
1134 faceToDualPoint_[facei] = meshMod.addPoint
1135 (
1136 mesh_.faceCentres()[facei],
1137 mesh_.faces()[facei][0], // masterPoint
1138 -1, // zoneID
1139 true // inCell
1140 );
1141 }
1142 // Detect whether different dual cells on either side of a face. This
1143 // would necessitate having a dual face built from the face and thus a
1144 // dual point at the face centre.
1145 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1146 {
1147 if (faceToDualPoint_[facei] == -1)
1148 {
1149 const face& f = mesh_.faces()[facei];
1150
1151 forAll(f, fp)
1152 {
1153 label ownDualCell = findDualCell(own[facei], f[fp]);
1154 label neiDualCell = findDualCell(nei[facei], f[fp]);
1155
1156 if (ownDualCell != neiDualCell)
1157 {
1158 faceToDualPoint_[facei] = meshMod.addPoint
1159 (
1160 mesh_.faceCentres()[facei],
1161 f[fp], // masterPoint
1162 -1, // zoneID
1163 true // inCell
1164 );
1165
1166 break;
1167 }
1168 }
1169 }
1170 }
1171
1172 // From edge to dual point
1173
1174 forAll(featureEdges, i)
1175 {
1176 label edgeI = featureEdges[i];
1177
1178 const edge& e = mesh_.edges()[edgeI];
1179
1180 edgeToDualPoint_[edgeI] = meshMod.addPoint
1181 (
1182 e.centre(mesh_.points()),
1183 e[0], // masterPoint
1184 -1, // zoneID
1185 true // inCell
1186 );
1187 }
1188
1189 // Detect whether different dual cells on either side of an edge. This
1190 // would neccesitate having a dual face built perpendicular to the edge
1191 // and thus a dual point at the mid of the edge.
1192 // Note: not really true - the face can be built without the edge centre!
1193 const labelListList& edgeCells = mesh_.edgeCells();
1194
1195 forAll(edgeCells, edgeI)
1196 {
1197 if (edgeToDualPoint_[edgeI] == -1)
1198 {
1199 const edge& e = mesh_.edges()[edgeI];
1200
1201 // We need a point on the edge if not all cells on both sides
1202 // are the same.
1203
1204 const labelList& eCells = mesh_.edgeCells()[edgeI];
1205
1206 label dualE0 = findDualCell(eCells[0], e[0]);
1207 label dualE1 = findDualCell(eCells[0], e[1]);
1208
1209 for (label i = 1; i < eCells.size(); i++)
1210 {
1211 label newDualE0 = findDualCell(eCells[i], e[0]);
1212
1213 if (dualE0 != newDualE0)
1214 {
1215 edgeToDualPoint_[edgeI] = meshMod.addPoint
1216 (
1217 e.centre(mesh_.points()),
1218 e[0], // masterPoint
1219 mesh_.pointZones().whichZone(e[0]), // zoneID
1220 true // inCell
1221 );
1222
1223 break;
1224 }
1225
1226 label newDualE1 = findDualCell(eCells[i], e[1]);
1227
1228 if (dualE1 != newDualE1)
1229 {
1230 edgeToDualPoint_[edgeI] = meshMod.addPoint
1231 (
1232 e.centre(mesh_.points()),
1233 e[1], // masterPoint
1234 mesh_.pointZones().whichZone(e[1]), // zoneID
1235 true // inCell
1236 );
1237
1238 break;
1239 }
1240 }
1241 }
1242 }
1243
1244 // Make sure all boundary edges emanating from feature points are
1245 // feature edges as well.
1246 forAll(singleCellFeaturePoints, i)
1247 {
1248 generateDualBoundaryEdges
1249 (
1250 isBoundaryEdge,
1251 singleCellFeaturePoints[i],
1252 meshMod
1253 );
1254 }
1255 forAll(multiCellFeaturePoints, i)
1256 {
1257 generateDualBoundaryEdges
1258 (
1259 isBoundaryEdge,
1260 multiCellFeaturePoints[i],
1261 meshMod
1262 );
1263 }
1264
1265
1266 // Check for duplicate points
1267 if (debug)
1268 {
1269 dumpPolyTopoChange(meshMod, "generatedPoints_");
1270 checkPolyTopoChange(meshMod);
1271 }
1272
1273
1274 // Now we have all points and cells
1275 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1276 // - pointToDualCells_ : per point a single dualCell or multiple dualCells
1277 // - pointToDualPoint_ : per point -1 or the dual point at the coordinate
1278 // - edgeToDualPoint_ : per edge -1 or the edge centre
1279 // - faceToDualPoint_ : per face -1 or the face centre
1280 // - cellToDualPoint_ : per cell the cell centre
1281 // Now we have to walk all edges and construct faces. Either single face
1282 // per edge or multiple (-if nonmanifold edge -if different dualcells)
1283
1284 const edgeList& edges = mesh_.edges();
1285
1286 forAll(edges, edgeI)
1287 {
1288 const labelList& eFaces = mesh_.edgeFaces()[edgeI];
1289
1290 boolList doneEFaces(eFaces.size(), false);
1291
1292 forAll(eFaces, i)
1293 {
1294 if (!doneEFaces[i])
1295 {
1296 // We found a face that hasn't yet been visited. This might
1297 // happen for non-manifold edges where a single edge can
1298 // become multiple faces.
1299
1300 label startFacei = eFaces[i];
1301
1302 //Pout<< "Walking edge:" << edgeI
1303 // << " points:" << mesh_.points()[e[0]]
1304 // << mesh_.points()[e[1]]
1305 // << " startFace:" << startFacei
1306 // << " at:" << mesh_.faceCentres()[startFacei]
1307 // << endl;
1308
1309 createFacesAroundEdge
1310 (
1311 splitFace,
1312 isBoundaryEdge,
1313 edgeI,
1314 startFacei,
1315 meshMod,
1316 doneEFaces
1317 );
1318 }
1319 }
1320 }
1321
1322 if (debug)
1323 {
1324 dumpPolyTopoChange(meshMod, "generatedFacesFromEdges_");
1325 }
1326
1327 // Create faces from feature faces. These can be internal or external faces.
1328 // - feature face : centre needs to be included.
1329 // - single cells on either side: triangulate
1330 // - multiple cells: create single face between unique cell pair. Only
1331 // create face where cells differ on either side.
1332 // - non-feature face : inbetween cell zones.
1333 forAll(faceToDualPoint_, facei)
1334 {
1335 if (faceToDualPoint_[facei] != -1 && mesh_.isInternalFace(facei))
1336 {
1337 const face& f = mesh_.faces()[facei];
1338 const labelList& fEdges = mesh_.faceEdges()[facei];
1339
1340 // Starting edge
1341 label fp = 0;
1342
1343 do
1344 {
1345 // Find edge that is in dual mesh and where the
1346 // next point (fp+1) has different dual cells on either side.
1347 bool foundStart = false;
1348
1349 do
1350 {
1351 if
1352 (
1353 edgeToDualPoint_[fEdges[fp]] != -1
1354 && !sameDualCell(facei, f.nextLabel(fp))
1355 )
1356 {
1357 foundStart = true;
1358 break;
1359 }
1360 fp = f.fcIndex(fp);
1361 }
1362 while (fp != 0);
1363
1364 if (!foundStart)
1365 {
1366 break;
1367 }
1368
1369 // Walk from edge fp and generate a face.
1370 createFaceFromInternalFace
1371 (
1372 facei,
1373 fp,
1374 meshMod
1375 );
1376 }
1377 while (fp != 0);
1378 }
1379 }
1380
1381 if (debug)
1382 {
1383 dumpPolyTopoChange(meshMod, "generatedFacesFromFeatFaces_");
1384 }
1385
1386
1387 // Create boundary faces. Every boundary point has one or more dualcells.
1388 // These need to be closed.
1389 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1390
1391 forAll(patches, patchi)
1392 {
1393 const polyPatch& pp = patches[patchi];
1394
1395 const labelListList& pointFaces = pp.pointFaces();
1396
1397 forAll(pointFaces, patchPointi)
1398 {
1399 const labelList& pFaces = pointFaces[patchPointi];
1400
1401 boolList donePFaces(pFaces.size(), false);
1402
1403 forAll(pFaces, i)
1404 {
1405 if (!donePFaces[i])
1406 {
1407 // Starting face
1408 label startFacei = pp.start()+pFaces[i];
1409
1410 //Pout<< "Walking around point:" << pointi
1411 // << " coord:" << mesh_.points()[pointi]
1412 // << " on patch:" << patchi
1413 // << " startFace:" << startFacei
1414 // << " at:" << mesh_.faceCentres()[startFacei]
1415 // << endl;
1416
1417 createFacesAroundBoundaryPoint
1418 (
1419 patchi,
1420 patchPointi,
1421 startFacei,
1422 meshMod,
1423 donePFaces // pFaces visited
1424 );
1425 }
1426 }
1427 }
1428 }
1429
1430 if (debug)
1431 {
1432 dumpPolyTopoChange(meshMod, "generatedFacesFromBndFaces_");
1433 }
1434}
1435
1436
1437// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
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
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const labelListList & pointFaces() const
Return point-face addressing.
A List with indirect addressing. Like IndirectList but does not store addressing.
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
Definition UListI.H:106
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
Walks from starting face around edge.
static label getMinIndex(const face &f, const label v0, const label v1)
Helper: find index in face of edge or -1. Index is such that edge is.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A 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
Creates dual of polyMesh. Every point becomes a cell (or multiple cells for feature points),...
void setRefinement(const bool splitFace, const labelList &featureFaces, const labelList &featureEdges, const labelList &singleCellFeaturePoints, const labelList &multiCellFeaturePoints, polyTopoChange &meshMod)
Insert all changes into meshMod to convert the polyMesh into.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Direct mesh changes based on v1.3 polyTopoChange syntax.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
Geometric merging of points. See below.
Namespace for handling debugging switches.
Definition debug.C:45
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
CompactListList< label > invertOneToManyCompact(const label len, const labelUList &map)
Invert one-to-many compact map. Unmapped elements will be size 0.
Definition ListOps.C:175
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
messageStream Info
Information stream (stdout output on master, null elsewhere).
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
Calculate merge mapping, preserving the original point order. All points closer/equal mergeTol are to...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
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.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299