conformalVoronoiMeshFeaturePoints.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-2015 OpenFOAM Foundation
9 Copyright (C) 2019 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 "vectorTools.H"
31#include "triangle.H"
32#include "tetrahedron.H"
33#include "Circulator.H"
34#include "DelaunayMeshTools.H"
35#include "OBJstream.H"
36
37using namespace Foam::vectorTools;
38
39// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40
42(
43 const extendedFeatureEdgeMesh& feMesh,
44 const pointIndexHit& edHit,
46) const
47{
48 if (foamyHexMeshControls().circulateEdges())
49 {
50 createEdgePointGroupByCirculating(feMesh, edHit, pts);
51 }
52 else
53 {
54 label edgeI = edHit.index();
55
57 feMesh.getEdgeStatus(edgeI);
58
59 switch (edStatus)
60 {
62 {
63 createExternalEdgePointGroup(feMesh, edHit, pts);
64 break;
65 }
67 {
68 createInternalEdgePointGroup(feMesh, edHit, pts);
69 break;
70 }
72 {
73 createFlatEdgePointGroup(feMesh, edHit, pts);
74 break;
75 }
77 {
78 createOpenEdgePointGroup(feMesh, edHit, pts);
79 break;
80 }
82 {
83 createMultipleEdgePointGroup(feMesh, edHit, pts);
84 break;
85 }
86 case extendedFeatureEdgeMesh::NONE:
87 {
88 break;
89 }
90 }
91 }
92}
93
94
95bool Foam::conformalVoronoiMesh::meshableRegion
96(
97 const plane::side side,
99) const
100{
101 switch (volType)
102 {
104 {
105 return (side == plane::BACK);
106 }
108 {
109 return (side == plane::FRONT);
110 }
112 {
113 return true;
114 }
116 {
117 return false;
118 }
119 default:
120 {
121 return false;
122 }
123 }
124}
125
126
127bool Foam::conformalVoronoiMesh::regionIsInside
128(
130 const vector& normalA,
132 const vector& normalB,
133 const vector& masterPtVec
134) const
135{
136 plane::side sideA
137 (
138 ((masterPtVec & normalA) <= 0) ? plane::BACK : plane::FRONT
139 );
140
141 plane::side sideB
142 (
143 ((masterPtVec & normalB) <= 0) ? plane::BACK : plane::FRONT
144 );
145
146 const bool meshableRegionA = meshableRegion(sideA, volTypeA);
147 const bool meshableRegionB = meshableRegion(sideB, volTypeB);
148
149 if (meshableRegionA == meshableRegionB)
150 {
151 return meshableRegionA;
152 }
153
155
156 return false;
157}
158
159
160void Foam::conformalVoronoiMesh::createEdgePointGroupByCirculating
161(
162 const extendedFeatureEdgeMesh& feMesh,
163 const pointIndexHit& edHit,
164 DynamicList<Vb>& pts
165) const
166{
167 typedef Foam::indexedVertexEnum::vertexType vertexType;
168 typedef extendedFeatureEdgeMesh::sideVolumeType sideVolumeType;
169
170 const Foam::point& edgePt = edHit.hitPoint();
171 const label edgeI = edHit.index();
172
173 scalar ppDist = pointPairDistance(edgePt);
174
175 const vectorField& feNormals = feMesh.normals();
176 const vector& edDir = feMesh.edgeDirections()[edgeI];
177 const labelList& edNormalIs = feMesh.edgeNormals()[edgeI];
178 const labelList& feNormalDirections = feMesh.normalDirections()[edgeI];
179
180 const List<sideVolumeType>& normalVolumeTypes = feMesh.normalVolumeTypes();
181
182 ConstCirculator<labelList> circ(edNormalIs);
183 ConstCirculator<labelList> circNormalDirs(feNormalDirections);
184
185 Map<Foam::point> masterPoints;
186 Map<vertexType> masterPointsTypes;
187 Map<plane> masterPointReflectionsPrev;
188 Map<plane> masterPointReflectionsNext;
189
190// Info<< "Edge = " << edHit << ", edDir = " << edDir << endl;
191// Info<< " edNorms = " << edNormalIs.size() << endl;
192
193 bool addedMasterPreviously = false;
194 label initialRegion = -1;
195
196 if (circ.size()) do
197 {
198 const sideVolumeType volType = normalVolumeTypes[circ.curr()];
199 const sideVolumeType nextVolType = normalVolumeTypes[circ.next()];
200
201 const vector& normal = feNormals[circ.curr()];
202 const vector& nextNormal = feNormals[circ.next()];
203
204 vector normalDir = (normal ^ edDir);
205 normalDir *= circNormalDirs.curr()/mag(normalDir);
206
207 vector nextNormalDir = (nextNormal ^ edDir);
208 nextNormalDir *= circNormalDirs.next()/mag(nextNormalDir);
209
210// Info<< " " << circ.curr() << " " << circ.next() << nl
211// << " " << circNormalDirs.curr() << " " << circNormalDirs.next()
212// << nl
213// << " normals = " << normal << " " << nextNormal << nl
214// << " normalDirs = " << normalDir << " " << nextNormalDir << nl
215// << " cross = " << (normalDir ^ nextNormalDir) << nl
216// << " "
217// << extendedFeatureEdgeMesh::sideVolumeTypeNames_[volType] << " "
218// << extendedFeatureEdgeMesh::sideVolumeTypeNames_[nextVolType]
219// << endl;
220
221 // Calculate master point
222 const vector masterPtVec = normalised(normalDir + nextNormalDir);
223
224 if
225 (
226 ((normalDir ^ nextNormalDir) & edDir) < SMALL
227 || mag(masterPtVec) < SMALL
228 )
229 {
230// Info<< " IGNORE REGION" << endl;
231 addedMasterPreviously = false;
232
233 if
234 (
235 circ.size() == 2
236 && mag((normal & nextNormal) - 1) < SMALL
237 )
238 {
239 const vector n = 0.5*(normal + nextNormal);
240
241 const vector s = ppDist*(edDir ^ n);
242
243 if (volType == extendedFeatureEdgeMesh::BOTH)
244 {
245 createBafflePointPair(ppDist, edgePt + s, n, true, pts);
246 createBafflePointPair(ppDist, edgePt - s, n, true, pts);
247 }
248 else
249 {
251 << "Failed to insert flat/open edge as volType is "
253 [
254 volType
255 ]
256 << endl;
257 }
258
259 break;
260 }
261
262 continue;
263 }
264
265 const Foam::point masterPt = edgePt + ppDist*masterPtVec;
266
267 // Check that region is inside or outside
268 const bool inside =
269 regionIsInside
270 (
271 volType,
272 normal,
273 nextVolType,
274 nextNormal,
275 masterPtVec
276 );
277
278 // Specialise for size = 1 && baffle
279 if (mag((normalDir & nextNormalDir) - 1) < SMALL)
280 {
281 if (inside)
282 {
283// Info<< "Specialise for size 1 and baffle" << endl;
284
285 vector s = ppDist*(edDir ^ normal);
286
287 plane facePlane(edgePt, normal);
288
289 Foam::point pt1 = edgePt + s + ppDist*normal;
290 Foam::point pt2 = edgePt - s + ppDist*normal;
291
292 Foam::point pt3 = facePlane.mirror(pt1);
293 Foam::point pt4 = facePlane.mirror(pt2);
294
299
300 break;
301 }
302 else
303 {
305 << "Faces are parallel but master point is not inside"
306 << endl;
307 }
308 }
309
310 if (!addedMasterPreviously)
311 {
312 if (initialRegion == -1)
313 {
314 initialRegion = circ.nRotations();
315 }
316
317 addedMasterPreviously = true;
318
319 masterPoints.insert(circ.nRotations(), masterPt);
320 masterPointsTypes.insert
321 (
322 circ.nRotations(),
323 inside
326 );
327
328 masterPointReflectionsPrev.insert
329 (
330 circ.nRotations(),
331 plane(edgePt, normal)
332 );
333
334 masterPointReflectionsNext.insert
335 (
336 circ.nRotations(),
337 plane(edgePt, nextNormal)
338 );
339 }
340 else if (addedMasterPreviously)
341 {
342 addedMasterPreviously = true;
343
344 masterPointReflectionsNext.erase(circ.nRotations() - 1);
345
346 // Shift the master point to be normal to the plane between it and
347 // the previous master point
348 // Should be the intersection of the normal and the plane with the
349 // new master point in it.
350
351 plane p(masterPoints[circ.nRotations() - 1], normalDir);
352 plane::ray r(edgePt, masterPt - edgePt);
353
354 scalar cutPoint = p.normalIntersect(r);
355
356 masterPoints.insert
357 (
358 circ.nRotations(),
359 edgePt + cutPoint*(masterPt - edgePt)
360 );
361
362 masterPointsTypes.insert
363 (
364 circ.nRotations(),
365 inside
368 );
369
370 masterPointReflectionsNext.insert
371 (
372 circ.nRotations(),
373 plane(edgePt, nextNormal)
374 );
375 }
376
377 if
378 (
379 masterPoints.size() > 1
380 && inside
381 && circ.nRotations() == circ.size() - 1
382 )
383 {
384 if (initialRegion == 0)
385 {
386 plane p(masterPoints[initialRegion], nextNormalDir);
387 plane::ray r(edgePt, masterPt - edgePt);
388
389 scalar cutPoint = p.normalIntersect(r);
390
391 masterPoints[circ.nRotations()] =
392 edgePt + cutPoint*(masterPt - edgePt);
393
394 // Remove the first reflection plane if we are no longer
395 // circulating
396
397 masterPointReflectionsPrev.erase(initialRegion);
398 masterPointReflectionsNext.erase(circ.nRotations());
399 }
400 else
401 {
402
403 }
404 }
405 }
406 while
407 (
408 circ.circulate(CirculatorBase::CLOCKWISE),
409 circNormalDirs.circulate(CirculatorBase::CLOCKWISE)
410 );
411
412
413 forAllConstIters(masterPoints, iter)
414 {
415 const Foam::point& pt = masterPoints[iter.key()];
416 const vertexType ptType = masterPointsTypes[iter.key()];
417
418// Info<< " Adding Master " << iter.key() << " " << pt << " "
419// << indexedVertexEnum::vertexTypeNames_[ptType] << endl;
420
421 pts.append(Vb(pt, ptType));
422
423 const vertexType reflectedPtType =
424 (
428 );
429
430 if (masterPointReflectionsPrev.found(iter.key()))
431 {
432 const Foam::point reflectedPt =
433 masterPointReflectionsPrev[iter.key()].mirror(pt);
434
435// Info<< " Adding Prev " << reflectedPt << " "
436// << indexedVertexEnum::vertexTypeNames_[reflectedPtType]
437// << endl;
438
439 pts.append(Vb(reflectedPt, reflectedPtType));
440 }
441
442 if (masterPointReflectionsNext.found(iter.key()))
443 {
444 const Foam::point reflectedPt =
445 masterPointReflectionsNext[iter.key()].mirror(pt);
446
447// Info<< " Adding Next " << reflectedPt << " "
448// << indexedVertexEnum::vertexTypeNames_[reflectedPtType]
449// << endl;
450
451 pts.append(Vb(reflectedPt, reflectedPtType));
452 }
453 }
454
455// pts.append(Vb(edgePt, Vb::vtExternalFeatureEdge));
456}
457
458
459void Foam::conformalVoronoiMesh::createExternalEdgePointGroup
460(
461 const extendedFeatureEdgeMesh& feMesh,
462 const pointIndexHit& edHit,
463 DynamicList<Vb>& pts
464) const
465{
466 const Foam::point& edgePt = edHit.hitPoint();
467
468 scalar ppDist = pointPairDistance(edgePt);
469
470 const vectorField& feNormals = feMesh.normals();
471 const labelList& edNormalIs = feMesh.edgeNormals()[edHit.index()];
472 const List<extendedFeatureEdgeMesh::sideVolumeType>& normalVolumeTypes =
473 feMesh.normalVolumeTypes();
474
475 // As this is an external edge, there are two normals by definition
476 const vector& nA = feNormals[edNormalIs[0]];
477 const vector& nB = feNormals[edNormalIs[1]];
478
480 normalVolumeTypes[edNormalIs[0]];
481
483 normalVolumeTypes[edNormalIs[1]];
484
485 if (areParallel(nA, nB))
486 {
487 // The normals are nearly parallel, so this is too sharp a feature to
488 // conform to.
489 return;
490 }
491
492 // Normalised distance of reference point from edge point
493 vector refVec((nA + nB)/(1 + (nA & nB)));
494
495 if (magSqr(refVec) > sqr(5.0))
496 {
497 // Limit the size of the conformation
498 ppDist *= 5.0/mag(refVec);
499
500 // Pout<< nl << "createExternalEdgePointGroup limit "
501 // << "edgePt " << edgePt << nl
502 // << "refVec " << refVec << nl
503 // << "mag(refVec) " << mag(refVec) << nl
504 // << "ppDist " << ppDist << nl
505 // << "nA " << nA << nl
506 // << "nB " << nB << nl
507 // << "(nA & nB) " << (nA & nB) << nl
508 // << endl;
509 }
510
511 // Convex. So refPt will be inside domain and hence a master point
512 Foam::point refPt = edgePt - ppDist*refVec;
513
514 // Insert the master point pairing the first slave
515
516 if (!geometryToConformTo_.inside(refPt))
517 {
518 return;
519 }
520
521 pts.append
522 (
523 Vb
524 (
525 refPt,
526 vertexCount() + pts.size(),
529 )
530 );
531
532 // Insert the slave points by reflecting refPt in both faces.
533 // with each slave referring to the master
534
535 Foam::point reflectedA = refPt + 2*ppDist*nA;
536 pts.append
537 (
538 Vb
539 (
540 reflectedA,
541 vertexCount() + pts.size(),
542 (
546 ),
548 )
549 );
550
551 Foam::point reflectedB = refPt + 2*ppDist*nB;
552 pts.append
553 (
554 Vb
555 (
556 reflectedB,
557 vertexCount() + pts.size(),
558 (
562 ),
564 )
565 );
566
567 ptPairs_.addPointPair
568 (
569 pts[pts.size() - 3].index(),
570 pts[pts.size() - 1].index()
571 );
572
573 ptPairs_.addPointPair
574 (
575 pts[pts.size() - 3].index(),
576 pts[pts.size() - 2].index()
577 );
578}
579
580
581void Foam::conformalVoronoiMesh::createInternalEdgePointGroup
582(
583 const extendedFeatureEdgeMesh& feMesh,
584 const pointIndexHit& edHit,
585 DynamicList<Vb>& pts
586) const
587{
588 const Foam::point& edgePt = edHit.hitPoint();
589
590 scalar ppDist = pointPairDistance(edgePt);
591
592 const vectorField& feNormals = feMesh.normals();
593 const labelList& edNormalIs = feMesh.edgeNormals()[edHit.index()];
594 const List<extendedFeatureEdgeMesh::sideVolumeType>& normalVolumeTypes =
595 feMesh.normalVolumeTypes();
596
597 // As this is an external edge, there are two normals by definition
598 const vector& nA = feNormals[edNormalIs[0]];
599 const vector& nB = feNormals[edNormalIs[1]];
600
602 normalVolumeTypes[edNormalIs[0]];
603
604// const extendedFeatureEdgeMesh::sideVolumeType& volTypeB =
605// normalVolumeTypes[edNormalIs[1]];
606
607 if (areParallel(nA, nB))
608 {
609 // The normals are nearly parallel, so this is too sharp a feature to
610 // conform to.
611 return;
612 }
613
614 // Normalised distance of reference point from edge point
615 vector refVec((nA + nB)/(1 + (nA & nB)));
616
617 if (magSqr(refVec) > sqr(5.0))
618 {
619 // Limit the size of the conformation
620 ppDist *= 5.0/mag(refVec);
621
622 // Pout<< nl << "createInternalEdgePointGroup limit "
623 // << "edgePt " << edgePt << nl
624 // << "refVec " << refVec << nl
625 // << "mag(refVec) " << mag(refVec) << nl
626 // << "ppDist " << ppDist << nl
627 // << "nA " << nA << nl
628 // << "nB " << nB << nl
629 // << "(nA & nB) " << (nA & nB) << nl
630 // << endl;
631 }
632
633 // Concave. master and reflected points inside the domain.
634 Foam::point refPt = edgePt - ppDist*refVec;
635
636 // Generate reflected master to be outside.
637 Foam::point reflMasterPt = refPt + 2*(edgePt - refPt);
638
639 // Reflect reflMasterPt in both faces.
640 Foam::point reflectedA = reflMasterPt - 2*ppDist*nA;
641
642 Foam::point reflectedB = reflMasterPt - 2*ppDist*nB;
643
644 scalar totalAngle =
646
647 // Number of quadrants the angle should be split into
648 int nQuads = int(totalAngle/foamyHexMeshControls().maxQuadAngle()) + 1;
649
650 // The number of additional master points needed to obtain the
651 // required number of quadrants.
652 int nAddPoints = min(max(nQuads - 2, 0), 2);
653
654 // Add number_of_vertices() at insertion of first vertex to all numbers:
655 // Result for nAddPoints 1 when the points are eventually inserted
656 // pt index type
657 // reflectedA 0 2
658 // reflectedB 1 2
659 // reflMasterPt 2 0
660
661 // Result for nAddPoints 1 when the points are eventually inserted
662 // pt index type
663 // reflectedA 0 3
664 // reflectedB 1 3
665 // refPt 2 3
666 // reflMasterPt 3 0
667
668 // Result for nAddPoints 2 when the points are eventually inserted
669 // pt index type
670 // reflectedA 0 4
671 // reflectedB 1 4
672 // reflectedAa 2 4
673 // reflectedBb 3 4
674 // reflMasterPt 4 0
675
676 if
677 (
678 !geometryToConformTo_.inside(reflectedA)
679 || !geometryToConformTo_.inside(reflectedB)
680 )
681 {
682 return;
683 }
684
685 // Master A is inside.
686 pts.append
687 (
688 Vb
689 (
690 reflectedA,
691 vertexCount() + pts.size(),
694 )
695 );
696
697 // Master B is inside.
698 pts.append
699 (
700 Vb
701 (
702 reflectedB,
703 vertexCount() + pts.size(),
706 )
707 );
708
709 // Slave is outside.
710 pts.append
711 (
712 Vb
713 (
714 reflMasterPt,
715 vertexCount() + pts.size(),
716 (
720 ),
722 )
723 );
724
725 ptPairs_.addPointPair
726 (
727 pts[pts.size() - 2].index(),
728 pts[pts.size() - 1].index()
729 );
730
731 ptPairs_.addPointPair
732 (
733 pts[pts.size() - 3].index(),
734 pts[pts.size() - 1].index()
735 );
736
737 if (nAddPoints == 1)
738 {
739 // One additional point is the reflection of the slave point,
740 // i.e. the original reference point
741 pts.append
742 (
743 Vb
744 (
745 refPt,
746 vertexCount() + pts.size(),
749 )
750 );
751 }
752 else if (nAddPoints == 2)
753 {
754 Foam::point reflectedAa = refPt + ppDist*nB;
755 pts.append
756 (
757 Vb
758 (
759 reflectedAa,
760 vertexCount() + pts.size(),
763 )
764 );
765
766 Foam::point reflectedBb = refPt + ppDist*nA;
767 pts.append
768 (
769 Vb
770 (
771 reflectedBb,
772 vertexCount() + pts.size(),
775 )
776 );
777 }
778}
779
780
781void Foam::conformalVoronoiMesh::createFlatEdgePointGroup
782(
783 const extendedFeatureEdgeMesh& feMesh,
784 const pointIndexHit& edHit,
785 DynamicList<Vb>& pts
786) const
787{
788 const Foam::point& edgePt = edHit.hitPoint();
789
790 const scalar ppDist = pointPairDistance(edgePt);
791
792 const vectorField& feNormals = feMesh.normals();
793 const labelList& edNormalIs = feMesh.edgeNormals()[edHit.index()];
794 const List<extendedFeatureEdgeMesh::sideVolumeType>& normalVolumeTypes =
795 feMesh.normalVolumeTypes();
796
797 // As this is a flat edge, there are two normals by definition
798 const vector& nA = feNormals[edNormalIs[0]];
799 const vector& nB = feNormals[edNormalIs[1]];
800
801 // Average normal to remove any bias to one side, although as this
802 // is a flat edge, the normals should be essentially the same
803 const vector n = 0.5*(nA + nB);
804
805 // Direction along the surface to the control point, sense of edge
806 // direction not important, as +s and -s can be used because this
807 // is a flat edge
808 vector s = ppDist*(feMesh.edgeDirections()[edHit.index()] ^ n);
809
810 if (normalVolumeTypes[edNormalIs[0]] == extendedFeatureEdgeMesh::OUTSIDE)
811 {
812 createPointPair(ppDist, edgePt + s, -n, true, pts);
813 createPointPair(ppDist, edgePt - s, -n, true, pts);
814 }
815 else if (normalVolumeTypes[edNormalIs[0]] == extendedFeatureEdgeMesh::BOTH)
816 {
817 createBafflePointPair(ppDist, edgePt + s, n, true, pts);
818 createBafflePointPair(ppDist, edgePt - s, n, true, pts);
819 }
820 else
821 {
822 createPointPair(ppDist, edgePt + s, n, true, pts);
823 createPointPair(ppDist, edgePt - s, n, true, pts);
824 }
825}
826
827
828void Foam::conformalVoronoiMesh::createOpenEdgePointGroup
829(
830 const extendedFeatureEdgeMesh& feMesh,
831 const pointIndexHit& edHit,
832 DynamicList<Vb>& pts
833) const
834{
835 // Assume it is a baffle and insert flat edge point pairs
836 const Foam::point& edgePt = edHit.hitPoint();
837
838 const scalar ppDist = pointPairDistance(edgePt);
839
840 const vectorField& feNormals = feMesh.normals();
841 const labelList& edNormalIs = feMesh.edgeNormals()[edHit.index()];
842
843 if (edNormalIs.size() == 1)
844 {
845// Info<< "Inserting open edge point group around " << edgePt << endl;
846// Info<< " ppDist = " << ppDist << nl
847// << " edNormals = " << edNormalIs
848// << endl;
849
850 const vector& n = feNormals[edNormalIs[0]];
851
852 vector s = ppDist*(feMesh.edgeDirections()[edHit.index()] ^ n);
853
854 plane facePlane(edgePt, n);
855
856 const label initialPtsSize = pts.size();
857
858 if
859 (
860 !geometryToConformTo_.inside(edgePt)
861 )
862 {
863 return;
864 }
865
866 createBafflePointPair(ppDist, edgePt - s, n, true, pts);
867 createBafflePointPair(ppDist, edgePt + s, n, false, pts);
868
869 for (label ptI = initialPtsSize; ptI < pts.size(); ++ptI)
870 {
871 pts[ptI].type() = Vb::vtInternalFeatureEdge;
872 }
873 }
874 else
875 {
876 Info<< "NOT INSERTING OPEN EDGE POINT GROUP WITH MORE THAN 1 "
877 << "EDGE NORMAL, NOT IMPLEMENTED" << endl;
878 }
879}
880
881
882void Foam::conformalVoronoiMesh::createMultipleEdgePointGroup
883(
884 const extendedFeatureEdgeMesh& feMesh,
885 const pointIndexHit& edHit,
886 DynamicList<Vb>& pts
887) const
888{
889// Info<< "NOT INSERTING MULTIPLE EDGE POINT GROUP, NOT IMPLEMENTED" << endl;
890
891 const Foam::point& edgePt = edHit.hitPoint();
892
893 const scalar ppDist = pointPairDistance(edgePt);
894
895 const vector edDir = feMesh.edgeDirections()[edHit.index()];
896
897 const vectorField& feNormals = feMesh.normals();
898 const labelList& edNormalIs = feMesh.edgeNormals()[edHit.index()];
899 const labelList& normalDirs = feMesh.normalDirections()[edHit.index()];
900
901 const List<extendedFeatureEdgeMesh::sideVolumeType>& normalVolumeTypes =
902 feMesh.normalVolumeTypes();
903
904 labelList nNormalTypes(4, Zero);
905
906 forAll(edNormalIs, edgeNormalI)
907 {
909 normalVolumeTypes[edNormalIs[edgeNormalI]];
910
911 nNormalTypes[sType]++;
912 }
913
914 if (nNormalTypes[extendedFeatureEdgeMesh::BOTH] == 4)
915 {
916 label masterEdgeNormalIndex = -1;
917
918 forAll(edNormalIs, edgeNormalI)
919 {
921 normalVolumeTypes[edNormalIs[edgeNormalI]];
922
924 {
925 masterEdgeNormalIndex = edgeNormalI;
926 break;
927 }
928 }
929
930 const vector& n = feNormals[edNormalIs[masterEdgeNormalIndex]];
931
932 label nDir = normalDirs[masterEdgeNormalIndex];
933
934 vector normalDir =
935 (feNormals[edNormalIs[masterEdgeNormalIndex]] ^ edDir);
936 normalDir *= nDir/mag(normalDir);
937
938 Foam::point pt1 = edgePt + ppDist*normalDir + ppDist*n;
939 Foam::point pt2 = edgePt + ppDist*normalDir - ppDist*n;
940
941 plane plane3(edgePt, normalDir);
942
943 Foam::point pt3 = plane3.mirror(pt1);
944 Foam::point pt4 = plane3.mirror(pt2);
945
946 pts.append
947 (
948 Vb
949 (
950 pt1,
951 vertexCount() + pts.size(),
954 )
955 );
956 pts.append
957 (
958 Vb
959 (
960 pt2,
961 vertexCount() + pts.size(),
964 )
965 );
966
967 ptPairs_.addPointPair
968 (
969 pts[pts.size() - 2].index(), // external 0 -> slave
970 pts[pts.size() - 1].index()
971 );
972
973 pts.append
974 (
975 Vb
976 (
977 pt3,
978 vertexCount() + pts.size(),
981 )
982 );
983
984 ptPairs_.addPointPair
985 (
986 pts[pts.size() - 3].index(), // external 0 -> slave
987 pts[pts.size() - 1].index()
988 );
989
990 pts.append
991 (
992 Vb
993 (
994 pt4,
995 vertexCount() + pts.size(),
998 )
999 );
1000
1001 ptPairs_.addPointPair
1002 (
1003 pts[pts.size() - 3].index(), // external 0 -> slave
1004 pts[pts.size() - 1].index()
1005 );
1006
1007 ptPairs_.addPointPair
1008 (
1009 pts[pts.size() - 2].index(), // external 0 -> slave
1010 pts[pts.size() - 1].index()
1011 );
1012 }
1013 else if
1014 (
1015 nNormalTypes[extendedFeatureEdgeMesh::BOTH] == 1
1016 && nNormalTypes[extendedFeatureEdgeMesh::INSIDE] == 2
1017 )
1018 {
1019 label masterEdgeNormalIndex = -1;
1020
1021 forAll(edNormalIs, edgeNormalI)
1022 {
1024 normalVolumeTypes[edNormalIs[edgeNormalI]];
1025
1026 if (sType == extendedFeatureEdgeMesh::BOTH)
1027 {
1028 masterEdgeNormalIndex = edgeNormalI;
1029 break;
1030 }
1031 }
1032
1033 const vector& n = feNormals[edNormalIs[masterEdgeNormalIndex]];
1034
1035 label nDir = normalDirs[masterEdgeNormalIndex];
1036
1037 vector normalDir =
1038 (feNormals[edNormalIs[masterEdgeNormalIndex]] ^ edDir);
1039 normalDir *= nDir/mag(normalDir);
1040
1041 const label nextNormalI =
1042 (masterEdgeNormalIndex + 1) % edNormalIs.size();
1043 if ((normalDir & feNormals[edNormalIs[nextNormalI]]) > 0)
1044 {
1045 normalDir *= -1;
1046 }
1047
1048 Foam::point pt1 = edgePt + ppDist*normalDir + ppDist*n;
1049 Foam::point pt2 = edgePt + ppDist*normalDir - ppDist*n;
1050
1051 plane plane3(edgePt, normalDir);
1052
1053 Foam::point pt3 = plane3.mirror(pt1);
1054 Foam::point pt4 = plane3.mirror(pt2);
1055
1056 pts.append
1057 (
1058 Vb
1059 (
1060 pt1,
1061 vertexCount() + pts.size(),
1064 )
1065 );
1066 pts.append
1067 (
1068 Vb
1069 (
1070 pt2,
1071 vertexCount() + pts.size(),
1074 )
1075 );
1076
1077 ptPairs_.addPointPair
1078 (
1079 pts[pts.size() - 2].index(), // external 0 -> slave
1080 pts[pts.size() - 1].index()
1081 );
1082
1083 pts.append
1084 (
1085 Vb
1086 (
1087 pt3,
1088 vertexCount() + pts.size(),
1091 )
1092 );
1093
1094 ptPairs_.addPointPair
1095 (
1096 pts[pts.size() - 3].index(), // external 0 -> slave
1097 pts[pts.size() - 1].index()
1098 );
1099
1100 pts.append
1101 (
1102 Vb
1103 (
1104 pt4,
1105 vertexCount() + pts.size(),
1108 )
1109 );
1110
1111 ptPairs_.addPointPair
1112 (
1113 pts[pts.size() - 3].index(), // external 0 -> slave
1114 pts[pts.size() - 1].index()
1115 );
1116 }
1117
1118
1119// // As this is a flat edge, there are two normals by definition
1120// const vector& nA = feNormals[edNormalIs[0]];
1121// const vector& nB = feNormals[edNormalIs[1]];
1122//
1123// // Average normal to remove any bias to one side, although as this
1124// // is a flat edge, the normals should be essentially the same
1125// const vector n = 0.5*(nA + nB);
1126//
1127// // Direction along the surface to the control point, sense of edge
1128// // direction not important, as +s and -s can be used because this
1129// // is a flat edge
1130// vector s = ppDist*(feMesh.edgeDirections()[edHit.index()] ^ n);
1131//
1132// createBafflePointPair(ppDist, edgePt + s, n, true, pts);
1133// createBafflePointPair(ppDist, edgePt - s, n, true, pts);
1134}
1135
1136
1137void Foam::conformalVoronoiMesh::insertFeaturePoints(bool distribute)
1138{
1139 Info<< nl
1140 << "Inserting feature points" << endl;
1141
1142 const label preFeaturePointSize(number_of_vertices());
1143
1144 if (Pstream::parRun() && distribute)
1145 {
1146 ftPtConformer_.distribute(decomposition());
1147 }
1148
1149 const List<Vb>& ftPtVertices = ftPtConformer_.featurePointVertices();
1150
1151 // Insert the created points directly as already distributed.
1152 Map<label> oldToNewIndices =
1153 this->DelaunayMesh<Delaunay>::insertPoints(ftPtVertices, true);
1154
1155 ftPtConformer_.reIndexPointPairs(oldToNewIndices);
1156
1157 label nFeatureVertices = number_of_vertices() - preFeaturePointSize;
1158 reduce(nFeatureVertices, sumOp<label>());
1159
1160 Info<< " Inserted " << nFeatureVertices << " feature vertices" << endl;
1161}
1162
1163
1164//Foam::scalar Foam::conformalVoronoiMesh::pyramidVolume
1165//(
1166// const Foam::point& apex,
1167// const Foam::point& a,
1168// const Foam::point& b,
1169// const Foam::point& c,
1170// const bool printInfo
1171//) const
1172//{
1173// triPointRef tri(a, b, c);
1174//
1175// tetPointRef tet(tri.a(), tri.b(), tri.c(), apex);
1176//
1177// scalar volume = tet.mag();
1178//
1186//
1187// if (printInfo)
1188// {
1189// Info<< "Calculating volume of pyramid..." << nl
1190// << " Apex : " << apex << nl
1191// << " Point a : " << a << nl
1192// << " Point b : " << b << nl
1193// << " Point c : " << c << nl
1194// << " Center : " << tri.centre() << nl
1195// << " Volume : " << volume << endl;
1196// }
1197//
1198// return volume;
1199//}
1200
1201
1202//void Foam::conformalVoronoiMesh::createPyramidMasterPoint
1203//(
1204// const Foam::point& apex,
1205// const vectorField& edgeDirections,
1206// Foam::point& masterPoint,
1207// vectorField& norms
1208//) const
1209//{
1210// pointField basePoints(edgeDirections.size() + 1);
1211//
1212// forAll(edgeDirections, eI)
1213// {
1214// basePoints[eI] = edgeDirections[eI] + apex;
1215// }
1216//
1217// basePoints[edgeDirections.size() + 1] = apex;
1218//
1219// face f(identity(edgeDirections.size()));
1220//
1221// pyramidPointFaceRef p(f, apex);
1222//
1223// const scalar ppDist = pointPairDistance(apex);
1224//
1225//
1226// vector unitDir = f.centre();
1227// unitDir /= mag(unitDir);
1228//
1229// masterPoint = apex + ppDist*unitDir;
1230//
1231// norms.setSize(edgeDirections.size());
1232//
1233// forAll(norms, nI)
1234// {
1235// norms[nI] =
1236// }
1237//}
1238
1239
1240//void Foam::conformalVoronoiMesh::createConvexConcaveFeaturePoints
1241//(
1242// DynamicList<Foam::point>& pts,
1243// DynamicList<label>& indices,
1244// DynamicList<label>& types
1245//)
1246//{
1247// const PtrList<extendedFeatureEdgeMesh>& feMeshes
1248// (
1249// geometryToConformTo_.features()
1250// );
1251//
1252// forAll(feMeshes, i)
1253// {
1254// const extendedFeatureEdgeMesh& feMesh(feMeshes[i]);
1255//
1256// for
1257// (
1258// label ptI = feMesh.convexStart();
1259// ptI < feMesh.mixedStart();
1260// ptI++
1261// )
1262// {
1263// const Foam::point& apex = feMesh.points()[ptI];
1264//
1265// if (!positionOnThisProc(apex))
1266// {
1267// continue;
1268// }
1269//
1270// const vectorField& featPtEdgeDirections
1271// = feMesh.featurePointEdgeDirections(ptI);
1272//
1273// Foam::point masterPoint;
1274// vectorField tetNorms;
1275//
1276// createPyramidMasterPoint
1277// (
1278// apex,
1279// featPtEdgeDirections,
1280// masterPoint,
1281// tetNorms
1282// );
1283//
1284//
1285//
1286// // Result when the points are eventually inserted (example n = 4)
1287// // Add number_of_vertices() at insertion of first vertex to all
1288// // numbers:
1289// // pt index type
1290// // internalPt 0 1
1291// // externalPt0 1 0
1292// // externalPt1 2 0
1293// // externalPt2 3 0
1294// // externalPt3 4 0
1295//
1296// // Result when the points are eventually inserted (example n = 5)
1297// // Add number_of_vertices() at insertion of first vertex to all
1298// // numbers:
1299// // pt index type
1300// // internalPt0 0 5
1301// // internalPt1 1 5
1302// // internalPt2 2 5
1303// // internalPt3 3 5
1304// // internalPt4 4 5
1305// // externalPt 5 4
1306//
1307// if (geometryToConformTo_.inside(masterPoint))
1308// {
1309//
1310// }
1311// else
1312// {
1313//
1314// }
1315//
1316// pts.append(masterPoint);
1317// indices.append(0);
1318// types.append(1);
1319//
1320// label internalPtIndex = -1;
1321//
1322// forAll(tetNorms, nI)
1323// {
1324// const vector& n = tetNorms[nI];
1325//
1326// Foam::point reflectedPoint
1327// = reflectPoint(featPt, masterPoint, n);
1328//
1329// pts.append(reflectedPoint);
1330// indices.append(0);
1331// types.append(internalPtIndex--);
1332// }
1333// }
1334// }
1335//}
CGAL::indexedVertex< K > Vb
label n
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:88
Like Foam::Circulator, but with const-access iterators.
Definition: Circulator.H:326
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:440
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
label index() const noexcept
Return the hit index.
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
void createEdgePointGroup(const extendedFeatureEdgeMesh &feMesh, const pointIndexHit &edHit, DynamicList< Vb > &pts) const
Call the appropriate function to conform to an edge.
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
const labelListList & normalDirections() const
const vectorField & edgeDirections() const
Return the edgeDirection vectors.
const List< sideVolumeType > & normalVolumeTypes() const
Return.
const labelListList & edgeNormals() const
Return the indices of the normals that are adjacent to the.
@ FLAT
Neither concave or convex, on a flat surface.
@ OPEN
Only connected to a single face.
@ MULTIPLE
Multiply connected (connected to more than two faces)
static const Enum< sideVolumeType > sideVolumeTypeNames_
const vectorField & normals() const
Return the normals of the surfaces adjacent to the feature edges.
sideVolumeType
Normals point to the outside.
@ NEITHER
not sure when this may be used
edgeStatus getEdgeStatus(label edgeI) const
Return the edgeStatus of a specified edge.
@ BOTH
limit by both minimum and maximum values
Definition: limitFields.H:189
A reference point and direction.
Definition: plane.H:110
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:95
side
Side of the plane.
Definition: plane.H:100
@ FRONT
The front (positive normal) side of the plane.
Definition: plane.H:101
@ BACK
The back (negative normal) side of the plane.
Definition: plane.H:102
int myProcNo() const noexcept
Return processor number.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
Collection of functions for testing relationships between two vectors.
Definition: vectorTools.H:50
bool areParallel(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Test if a and b are parallel: a^b = 0.
Definition: vectorTools.H:56
T radAngleBetween(const Vector< T > &a, const Vector< T > &b, const T &tolerance=SMALL)
Calculate angle between a and b in radians.
Definition: vectorTools.H:122
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278