Loading...
Searching...
No Matches
snappySnapDriverFeature.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-2015 OpenFOAM Foundation
9 Copyright (C) 2015-2022,2024 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 "snappySnapDriver.H"
30#include "polyTopoChange.H"
31#include "syncTools.H"
32#include "fvMesh.H"
33#include "OBJstream.H"
34#include "motionSmoother.H"
35#include "refinementSurfaces.H"
36#include "refinementFeatures.H"
37#include "unitConversion.H"
38#include "plane.H"
39#include "featureEdgeMesh.H"
40#include "treeDataPoint.H"
41#include "indexedOctree.H"
42#include "snapParameters.H"
43#include "PatchTools.H"
44#include "pyramid.H"
46
47// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48
49namespace Foam
50{
51 template<class T>
52 class listPlusEqOp
53 {
54 public:
55
56 void operator()(List<T>& x, const List<T>& y) const
57 {
58 label sz = x.size();
59 x.setSize(sz+y.size());
60 forAll(y, i)
61 {
62 x[sz++] = y[i];
63 }
64 }
65 };
66}
67
68
69// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
70
71//bool Foam::snappySnapDriver::isFeaturePoint
72//(
73// const scalar featureCos,
74// const indirectPrimitivePatch& pp,
75// const bitSet& isFeatureEdge,
76// const label pointi
77//) const
78//{
79// const pointField& points = ppLocalPoints;
80// const edgeList& edges = pp.edges();
81// const labelList& pEdges = pp.pointEdges()[pointi];
82//
83// label nFeatEdges = 0;
84//
85// forAll(pEdges, i)
86// {
87// if (isFeatureEdge[pEdges[i]])
88// {
89// nFeatEdges++;
90//
91// for (label j = i+1; j < pEdges.size(); j++)
92// {
93// if (isFeatureEdge[pEdges[j]])
94// {
95// const edge& ei = edges[pEdges[i]];
96// const edge& ej = edges[pEdges[j]];
97//
98// const point& p = points[pointi];
99// const point& pi = points[ei.otherVertex(pointi)];
100// const point& pj = points[ej.otherVertex(pointi)];
101//
102// vector vi = p-pi;
103// scalar viMag = mag(vi);
104//
105// vector vj = pj-p;
106// scalar vjMag = mag(vj);
107//
108// if
109// (
110// viMag > SMALL
111// && vjMag > SMALL
112// && ((vi/viMag & vj/vjMag) < featureCos)
113// )
114// {
115// return true;
116// }
117// }
118// }
119// }
120// }
121//
122// if (nFeatEdges == 1)
123// {
124// // End of feature-edge string
125// return true;
126// }
127//
128// return false;
129//}
130
131
132void Foam::snappySnapDriver::smoothAndConstrain
133(
134 const bitSet& isPatchMasterEdge,
136 const labelList& meshEdges,
137 const List<pointConstraint>& constraints,
138 vectorField& disp
139) const
140{
141 const fvMesh& mesh = meshRefiner_.mesh();
142
143 for (label avgIter = 0; avgIter < 20; avgIter++)
144 {
145 // Calculate average displacement of neighbours
146 // - unconstrained (i.e. surface) points use average of all
147 // neighbouring points
148 // - from testing it has been observed that it is not beneficial
149 // to have edge constrained points use average of all edge or point
150 // constrained neighbours since they're already attracted to
151 // the nearest point on the feature.
152 // Having them attract to point-constrained neighbours does not
153 // make sense either since there is usually just one of them so
154 // it severely distorts it.
155 // - same for feature points. They are already attracted to the
156 // nearest feature point.
157
158 vectorField dispSum(pp.nPoints(), Zero);
159 labelList dispCount(pp.nPoints(), Zero);
160
161 const labelListList& pointEdges = pp.pointEdges();
162 const edgeList& edges = pp.edges();
163
164 forAll(pointEdges, pointi)
165 {
166 const labelList& pEdges = pointEdges[pointi];
167
168 label nConstraints = constraints[pointi].first();
169
170 if (nConstraints <= 1)
171 {
172 forAll(pEdges, i)
173 {
174 label edgei = pEdges[i];
175
176 if (isPatchMasterEdge[edgei])
177 {
178 label nbrPointi = edges[edgei].otherVertex(pointi);
179 if (constraints[nbrPointi].first() >= nConstraints)
180 {
181 dispSum[pointi] += disp[nbrPointi];
182 dispCount[pointi]++;
183 }
184 }
185 }
186 }
187 }
188
190 (
191 mesh,
192 pp.meshPoints(),
193 dispSum,
194 plusEqOp<point>(),
196 mapDistribute::transform()
197 );
199 (
200 mesh,
201 pp.meshPoints(),
202 dispCount,
203 plusEqOp<label>(),
204 label(0),
205 mapDistribute::transform()
206 );
207
208 // Constraints
209 forAll(constraints, pointi)
210 {
211 if (dispCount[pointi] > 0)
212 {
213 // Mix my displacement with neighbours' displacement
214 disp[pointi] =
215 0.5
216 *(disp[pointi] + dispSum[pointi]/dispCount[pointi]);
217 }
218 }
219 }
220}
221
222
223void Foam::snappySnapDriver::calcNearestFace
224(
225 const label iter,
227 const pointField& ppLocalPoints,
228
229 const scalarField& faceSnapDist,
230 vectorField& faceDisp,
231 vectorField& faceSurfaceNormal,
232 labelList& faceSurfaceGlobalRegion
233 //vectorField& faceRotation
234) const
235{
236 const fvMesh& mesh = meshRefiner_.mesh();
237 const refinementSurfaces& surfaces = meshRefiner_.surfaces();
238
239 const pointField ppFaceCentres
240 (
242 (
244 ppLocalPoints
245 ).faceCentres()
246 );
247
248
249 // Displacement and orientation per pp face.
250 faceDisp.setSize(pp.size());
251 faceDisp = Zero;
252 faceSurfaceNormal.setSize(pp.size());
253 faceSurfaceNormal = Zero;
254 faceSurfaceGlobalRegion.setSize(pp.size());
255 faceSurfaceGlobalRegion = -1;
256
257 // Divide surfaces into zoned and unzoned
258 const labelList zonedSurfaces =
259 surfaceZonesInfo::getNamedSurfaces(surfaces.surfZones());
260 const labelList unzonedSurfaces =
261 surfaceZonesInfo::getUnnamedSurfaces(surfaces.surfZones());
262
263 // Per pp face the current surface snapped to
264 labelList snapSurf(pp.size(), -1);
265
266
267 // Do zoned surfaces
268 // ~~~~~~~~~~~~~~~~~
269 // Zoned faces only attract to corresponding surface
270
271 // Extract faces per zone
272 const PtrList<surfaceZonesInfo>& surfZones = surfaces.surfZones();
273
274 forAll(zonedSurfaces, i)
275 {
276 label zoneSurfi = zonedSurfaces[i];
277
278 const wordList& faceZoneNames = surfZones[zoneSurfi].faceZoneNames();
279
280 // Get indices of faces on pp that are also in zone
281 DynamicList<label> ppFaces;
282 DynamicList<label> meshFaces;
284 forAll(faceZoneNames, fzi)
285 {
286 const word& faceZoneName = faceZoneNames[fzi];
287 label zonei = mesh.faceZones().findZoneID(faceZoneName);
288 if (zonei == -1)
289 {
291 << "Problem. Cannot find zone " << faceZoneName
292 << exit(FatalError);
293 }
294 const faceZone& fZone = mesh.faceZones()[zonei];
295 const bitSet isZonedFace(mesh.nFaces(), fZone);
296
297 ppFaces.reserve(ppFaces.capacity()+fZone.size());
298 meshFaces.reserve(meshFaces.capacity()+fZone.size());
299 fc.reserve(meshFaces.capacity()+fZone.size());
300
301 forAll(pp.addressing(), i)
302 {
303 if (isZonedFace[pp.addressing()[i]])
304 {
305 snapSurf[i] = zoneSurfi;
306 ppFaces.append(i);
307 meshFaces.append(pp.addressing()[i]);
308 fc.append(ppFaceCentres[i]);
309 }
310 }
311
312 //Pout<< "For faceZone " << fZone.name()
313 // << " found " << ppFaces.size() << " out of " << pp.size()
314 // << endl;
315 }
316
317 List<pointIndexHit> hitInfo;
318 labelList hitSurface;
319 labelList hitRegion;
320 vectorField hitNormal;
321 surfaces.findNearestRegion
322 (
323 labelList(1, zoneSurfi),
324 fc,
325 sqr(scalarField(faceSnapDist, ppFaces)),// sqr of attract dist
326 hitSurface,
327 hitInfo,
328 hitRegion,
329 hitNormal
330 );
331
332 forAll(hitInfo, hiti)
333 {
334 if (hitInfo[hiti].hit())
335 {
336 label facei = ppFaces[hiti];
337 faceDisp[facei] = hitInfo[hiti].point() - fc[hiti];
338 faceSurfaceNormal[facei] = hitNormal[hiti];
339 faceSurfaceGlobalRegion[facei] = surfaces.globalRegion
340 (
341 hitSurface[hiti],
342 hitRegion[hiti]
343 );
344 }
345 else
346 {
347 static label nWarn = 0;
348
349 if (nWarn < 100)
350 {
352 << "Did not find surface near face centre " << fc[hiti]
353 << endl;
354 nWarn++;
355 if (nWarn == 100)
356 {
358 << "Reached warning limit " << nWarn
359 << ". Suppressing further warnings." << endl;
360 }
361 }
362 }
363 }
364 }
365
366
367 // Do unzoned surfaces
368 // ~~~~~~~~~~~~~~~~~~~
369 // Unzoned faces attract to any unzoned surface
370
371 DynamicList<label> ppFaces(pp.size());
372 DynamicList<label> meshFaces(pp.size());
374 forAll(pp.addressing(), i)
375 {
376 if (snapSurf[i] == -1)
377 {
378 ppFaces.append(i);
379 meshFaces.append(pp.addressing()[i]);
380 fc.append(ppFaceCentres[i]);
381 }
382 }
383 //Pout<< "Found " << ppFaces.size() << " unzoned faces out of "
384 // << pp.size() << endl;
385
386 List<pointIndexHit> hitInfo;
387 labelList hitSurface;
388 labelList hitRegion;
389 vectorField hitNormal;
390 surfaces.findNearestRegion
391 (
392 unzonedSurfaces,
393 fc,
394 sqr(scalarField(faceSnapDist, ppFaces)),// sqr of attract dist
395 hitSurface,
396 hitInfo,
397 hitRegion,
398 hitNormal
399 );
400
401 forAll(hitInfo, hiti)
402 {
403 if (hitInfo[hiti].hit())
404 {
405 label facei = ppFaces[hiti];
406 faceDisp[facei] = hitInfo[hiti].point() - fc[hiti];
407 faceSurfaceNormal[facei] = hitNormal[hiti];
408 faceSurfaceGlobalRegion[facei] = surfaces.globalRegion
409 (
410 hitSurface[hiti],
411 hitRegion[hiti]
412 );
413 }
414 else
415 {
416 static label nWarn = 0;
417
418 if (nWarn < 100)
419 {
421 << "Did not find surface near face centre " << fc[hiti]
422 << endl;
423
424 nWarn++;
425 if (nWarn == 100)
426 {
428 << "Reached warning limit " << nWarn
429 << ". Suppressing further warnings." << endl;
430 }
431 }
432 }
433 }
434
435
438 //
440 //faceRotation.setSize(pp.size());
441 //faceRotation = Zero;
442 //
443 //forAll(faceRotation, facei)
444 //{
445 // // Note: extend to >180 degrees checking
446 // faceRotation[facei] =
447 // pp.faceNormals()[facei]
448 // ^ faceSurfaceNormal[facei];
449 //}
450 //
451 //if (debug&meshRefinement::ATTRACTION)
452 //{
453 // dumpMove
454 // (
455 // mesh.time().path()
456 // / "faceDisp_" + name(iter) + ".obj",
457 // pp.faceCentres(),
458 // pp.faceCentres() + faceDisp
459 // );
460 // dumpMove
461 // (
462 // mesh.time().path()
463 // / "faceRotation_" + name(iter) + ".obj",
464 // pp.faceCentres(),
465 // pp.faceCentres() + faceRotation
466 // );
467 //}
468}
469
470
471// Collect (possibly remote) per point data of all surrounding faces
472// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
473// - faceSurfaceNormal
474// - faceDisp
475// - faceCentres&faceNormal
476void Foam::snappySnapDriver::calcNearestFacePointProperties
477(
478 const label iter,
480 const pointField& ppLocalPoints,
481
482 const vectorField& faceDisp,
483 const vectorField& faceSurfaceNormal,
484 const labelList& faceSurfaceGlobalRegion,
485
486 List<List<point>>& pointFaceSurfNormals,
487 List<List<point>>& pointFaceDisp,
488 List<List<point>>& pointFaceCentres,
489 List<labelList>& pointFacePatchID
490) const
491{
492 const fvMesh& mesh = meshRefiner_.mesh();
493
494 const bitSet isMasterFace(syncTools::getMasterFaces(mesh));
495
496
497 // For now just get all surrounding face data. Expensive - should just
498 // store and sync data on coupled points only
499 // (see e.g PatchToolsNormals.C)
500
501 pointFaceSurfNormals.setSize(pp.nPoints());
502 pointFaceDisp.setSize(pp.nPoints());
503 pointFaceCentres.setSize(pp.nPoints());
504 pointFacePatchID.setSize(pp.nPoints());
505
506 // Fill local data
507 forAll(pp.pointFaces(), pointi)
508 {
509 const labelList& pFaces = pp.pointFaces()[pointi];
510
511 // Count valid face normals
512 label nFaces = 0;
513 forAll(pFaces, i)
514 {
515 label facei = pFaces[i];
516 label globalRegioni = faceSurfaceGlobalRegion[facei];
517
518 if (isMasterFace[pp.addressing()[facei]] && globalRegioni != -1)
519 {
520 nFaces++;
521 }
522 }
523
524
525 List<point>& pNormals = pointFaceSurfNormals[pointi];
526 pNormals.setSize(nFaces);
527 List<point>& pDisp = pointFaceDisp[pointi];
528 pDisp.setSize(nFaces);
529 List<point>& pFc = pointFaceCentres[pointi];
530 pFc.setSize(nFaces);
531 labelList& pFid = pointFacePatchID[pointi];
532 pFid.setSize(nFaces);
533
534 nFaces = 0;
535 forAll(pFaces, i)
536 {
537 label facei = pFaces[i];
538 label globalRegioni = faceSurfaceGlobalRegion[facei];
539
540 if (isMasterFace[pp.addressing()[facei]] && globalRegioni != -1)
541 {
542 pNormals[nFaces] = faceSurfaceNormal[facei];
543 pDisp[nFaces] = faceDisp[facei];
544 pFc[nFaces] = pp.localFaces()[facei].centre(ppLocalPoints);
545 pFid[nFaces] = globalToMasterPatch_[globalRegioni];
546 nFaces++;
547 }
548 }
549 }
550
551
552 // Collect additionally 'normal' boundary faces for boundaryPoints of pp
553 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
554 // points on the boundary of pp should pick up non-pp normals
555 // as well for the feature-reconstruction to behave correctly.
556 // (the movement is already constrained outside correctly so it
557 // is only that the unconstrained attraction vector is calculated
558 // correctly)
559 {
562
563 // Unmark all non-coupled boundary faces
564 forAll(pbm, patchi)
565 {
566 const polyPatch& pp = pbm[patchi];
567
568 if (pp.coupled() || isA<emptyPolyPatch>(pp))
569 {
570 forAll(pp, i)
571 {
572 label meshFacei = pp.start()+i;
573 patchID[meshFacei-mesh.nInternalFaces()] = -1;
574 }
575 }
576 }
577
578 // Remove any meshed faces
579 forAll(pp.addressing(), i)
580 {
581 label meshFacei = pp.addressing()[i];
582 patchID[meshFacei-mesh.nInternalFaces()] = -1;
583 }
584
585
586
587 // See if edge of pp uses any non-meshed boundary faces. If so add the
588 // boundary face as additional constraint. Note that we account for
589 // both 'real' boundary edges and boundary edge of baffles
590
591 const labelList bafflePair
592 (
594 );
595
596
597 // Mark all points on 'boundary' edges
598 bitSet isBoundaryPoint(pp.nPoints());
599
600 const labelListList& edgeFaces = pp.edgeFaces();
601 const edgeList& edges = pp.edges();
602
603 forAll(edgeFaces, edgei)
604 {
605 const edge& e = edges[edgei];
606 const labelList& eFaces = edgeFaces[edgei];
607
608 if (eFaces.size() == 1)
609 {
610 // 'real' boundary edge
611 isBoundaryPoint.set(e[0]);
612 isBoundaryPoint.set(e[1]);
613 }
614 else if (eFaces.size() == 2 && bafflePair[eFaces[0]] == eFaces[1])
615 {
616 // 'baffle' boundary edge
617 isBoundaryPoint.set(e[0]);
618 isBoundaryPoint.set(e[1]);
619 }
620 }
621
622
623 // Construct labelList equivalent of meshPointMap
624 labelList meshToPatchPoint(mesh.nPoints(), -1);
625 forAll(pp.meshPoints(), pointi)
626 {
627 meshToPatchPoint[pp.meshPoints()[pointi]] = pointi;
628 }
629
630 forAll(patchID, bFacei)
631 {
632 label patchi = patchID[bFacei];
633
634 if (patchi != -1)
635 {
636 label facei = mesh.nInternalFaces()+bFacei;
637 const face& f = mesh.faces()[facei];
638
639 forAll(f, fp)
640 {
641 label pointi = meshToPatchPoint[f[fp]];
642
643 if (pointi != -1 && isBoundaryPoint.test(pointi))
644 {
645 List<point>& pNormals = pointFaceSurfNormals[pointi];
646 List<point>& pDisp = pointFaceDisp[pointi];
647 List<point>& pFc = pointFaceCentres[pointi];
648 labelList& pFid = pointFacePatchID[pointi];
649
650
651 // Note: these are points which are on the patch but
652 // the face itself is not. Could recalculate the
653 // facearea with the modified point but hopefully
654 // not necessary.
655 //const point& pt = mesh.points()[f[fp]];
656 const point& pt = ppLocalPoints[pointi];
657 vector fn = mesh.faceAreas()[facei];
658
659 pNormals.append(fn/mag(fn));
660 pDisp.append(mesh.faceCentres()[facei]-pt);
661 pFc.append(mesh.faceCentres()[facei]);
662 pFid.append(patchi);
663 }
664 }
665 }
666 }
667 }
668
670 (
671 mesh,
672 pp.meshPoints(),
673 pointFaceSurfNormals,
675 List<point>(),
677 );
679 (
680 mesh,
681 pp.meshPoints(),
682 pointFaceDisp,
684 List<point>(),
686 );
687
688 {
689 // Make into displacement before synchronising to avoid any problems
690 // with parallel cyclics
691 //pointField localPoints(pp.points(), pp.meshPoints());
692 forAll(pointFaceCentres, pointi)
693 {
694 //const point& pt = pp.points()[pp.meshPoints()[pointi]];
695 const point& pt = ppLocalPoints[pointi];
696
697 List<point>& pFc = pointFaceCentres[pointi];
698 for (point& p : pFc)
699 {
700 p -= pt;
701 }
702 }
704 (
705 mesh,
706 pp.meshPoints(),
707 pointFaceCentres,
709 List<point>(),
711 );
712 forAll(pointFaceCentres, pointi)
713 {
714 //const point& pt = pp.points()[pp.meshPoints()[pointi]];
715 const point& pt = ppLocalPoints[pointi];
716
717 List<point>& pFc = pointFaceCentres[pointi];
718 for (point& p : pFc)
719 {
720 p += pt;
721 }
722 }
723 }
724
726 (
727 mesh,
728 pp.meshPoints(),
729 pointFacePatchID,
732 );
733
734
735 // Sort the data according to the face centres. This is only so we get
736 // consistent behaviour serial and parallel.
737 labelList visitOrder;
738 forAll(pointFaceDisp, pointi)
739 {
740 List<point>& pNormals = pointFaceSurfNormals[pointi];
741 List<point>& pDisp = pointFaceDisp[pointi];
742 List<point>& pFc = pointFaceCentres[pointi];
743 labelList& pFid = pointFacePatchID[pointi];
744
745 sortedOrder(mag(pFc)(), visitOrder);
746
747 pNormals = List<point>(pNormals, visitOrder);
748 pDisp = List<point>(pDisp, visitOrder);
749 pFc = List<point>(pFc, visitOrder);
750 pFid = labelUIndList(pFid, visitOrder)();
751 }
752}
753
754
755// Gets passed in offset to nearest point on feature edge. Calculates
756// if the point has a different number of faces on either side of the feature
757// and if so attracts the point to that non-dominant plane.
758void Foam::snappySnapDriver::correctAttraction
759(
760 const UList<point>& surfacePoints,
761 const UList<label>& surfaceCounts,
762 const point& edgePt,
763 const vector& edgeNormal, // normalised normal
764 const point& pt,
765
766 vector& edgeOffset // offset from pt to point on edge
767) const
768{
769 // Tangential component along edge
770 scalar tang = ((pt-edgePt)&edgeNormal);
771
772 labelList order(sortedOrder(surfaceCounts));
773
774 if (order[0] < order[1])
775 {
776 // There is a non-dominant plane. Use the point on the plane to
777 // attract to.
778 vector attractD = surfacePoints[order[0]]-edgePt;
779 // Tangential component along edge
780 scalar tang2 = (attractD&edgeNormal);
781 // Normal component
782 attractD -= tang2*edgeNormal;
783 // Calculate fraction of normal distances
784 scalar magAttractD = mag(attractD);
785 scalar fraction = magAttractD/(magAttractD+mag(edgeOffset));
786
787 point linePt =
788 edgePt
789 + ((1.0-fraction)*tang2 + fraction*tang)*edgeNormal;
790 edgeOffset = linePt-pt;
791 }
792}
793
794
795Foam::pointIndexHit Foam::snappySnapDriver::findMultiPatchPoint
796(
797 const point& pt,
798 const labelList& patchIDs,
799 const List<point>& faceCentres
800) const
801{
802 // Determine if multiple patchIDs
803 if (patchIDs.size())
804 {
805 label patch0 = patchIDs[0];
806
807 for (label i = 1; i < patchIDs.size(); i++)
808 {
809 if (patchIDs[i] != patch0)
810 {
811 return pointIndexHit(true, pt, labelMax);
812 }
813 }
814 }
815 return pointIndexHit(false, Zero, labelMax);
816}
817
818
819Foam::label Foam::snappySnapDriver::findNormal
820(
821 const scalar featureCos,
822 const vector& n,
823 const UList<vector>& surfaceNormals
824) const
825{
826 label index = -1;
827
828 forAll(surfaceNormals, j)
829 {
830 scalar cosAngle = (n&surfaceNormals[j]);
831
832 if
833 (
834 (cosAngle >= featureCos)
835 || (cosAngle < (-1+0.001)) // triangle baffles
836 )
837 {
838 index = j;
839 break;
840 }
841 }
842 return index;
843}
844
845
846// Detect multiple patches. Returns pointIndexHit:
847// - false, index=-1 : single patch
848// - true , index=0 : multiple patches but on different normals planes
849// (so geometric feature edge is also a region edge)
850// - true , index=1 : multiple patches on same normals plane i.e. flat region
851// edge
852Foam::pointIndexHit Foam::snappySnapDriver::findMultiPatchPoint
853(
854 const point& pt,
855 const labelList& patchIDs,
856 const UList<vector>& surfaceNormals,
857 const labelList& faceToNormalBin
858) const
859{
860 if (patchIDs.empty())
861 {
862 return pointIndexHit(false, pt, -1);
863 }
864
865 // Detect single patch situation (to avoid allocation)
866 label patch0 = patchIDs[0];
867
868 for (label i = 1; i < patchIDs.size(); i++)
869 {
870 if (patchIDs[i] != patch0)
871 {
872 patch0 = -1;
873 break;
874 }
875 }
876
877 if (patch0 >= 0)
878 {
879 // Single patch
880 return pointIndexHit(false, pt, -1);
881 }
882 else
883 {
884 if (surfaceNormals.size() == 1)
885 {
886 // Same normals plane, flat region edge.
887 return pointIndexHit(true, pt, 1);
888 }
889 else
890 {
891 // Detect per normals bin
892 labelList normalToPatch(surfaceNormals.size(), -1);
893 forAll(faceToNormalBin, i)
894 {
895 if (faceToNormalBin[i] != -1)
896 {
897 label& patch = normalToPatch[faceToNormalBin[i]];
898 if (patch == -1)
899 {
900 // First occurrence
901 patch = patchIDs[i];
902 }
903 else if (patch == -2)
904 {
905 // Already marked as being on multiple patches
906 }
907 else if (patch != patchIDs[i])
908 {
909 // Mark as being on multiple patches
910 patch = -2;
911 }
912 }
913 }
914
915 forAll(normalToPatch, normali)
916 {
917 if (normalToPatch[normali] == -2)
918 {
919 // Multiple patches on same normals plane, flat region
920 // edge
921 return pointIndexHit(true, pt, 1);
922 }
923 }
924
925 // All patches on either side of geometric feature anyway
926 return pointIndexHit(true, pt, 0);
927 }
928 }
929}
930
931
932void Foam::snappySnapDriver::writeStats
933(
935 const bitSet& isPatchMasterPoint,
936 const List<pointConstraint>& patchConstraints
937) const
938{
939 label nMasterPoints = 0;
940 label nPlanar = 0;
941 label nEdge = 0;
942 label nPoint = 0;
943
944 forAll(patchConstraints, pointi)
945 {
946 if (isPatchMasterPoint[pointi])
947 {
948 nMasterPoints++;
949
950 if (patchConstraints[pointi].first() == 1)
951 {
952 nPlanar++;
953 }
954 else if (patchConstraints[pointi].first() == 2)
955 {
956 nEdge++;
957 }
958 else if (patchConstraints[pointi].first() == 3)
959 {
960 nPoint++;
961 }
962 }
963 }
964
965 reduce(nMasterPoints, sumOp<label>());
966 reduce(nPlanar, sumOp<label>());
967 reduce(nEdge, sumOp<label>());
968 reduce(nPoint, sumOp<label>());
969 Info<< "total master points :" << nMasterPoints
970 << " of which attracted to :" << nl
971 << " feature point : " << nPoint << nl
972 << " feature edge : " << nEdge << nl
973 << " nearest surface : " << nPlanar << nl
974 << " rest : " << nMasterPoints-nPoint-nEdge-nPlanar
975 << nl
976 << endl;
977}
978
979
980void Foam::snappySnapDriver::featureAttractionUsingReconstruction
981(
982 const label iter,
983 const scalar featureCos,
984
986 const pointField& ppLocalPoints,
987 const scalarField& snapDist,
988 const vectorField& nearestDisp,
989 const label pointi,
990
991 const List<List<point>>& pointFaceSurfNormals,
992 const List<List<point>>& pointFaceDisp,
993 const List<List<point>>& pointFaceCentres,
994 const labelListList& pointFacePatchID,
995
996 DynamicList<point>& surfacePoints,
997 DynamicList<vector>& surfaceNormals,
998 labelList& faceToNormalBin,
999
1000 vector& patchAttraction,
1001 pointConstraint& patchConstraint
1002) const
1003{
1004 patchAttraction = Zero;
1005 patchConstraint = pointConstraint();
1006
1007 const List<point>& pfSurfNormals = pointFaceSurfNormals[pointi];
1008 const List<point>& pfDisp = pointFaceDisp[pointi];
1009 const List<point>& pfCentres = pointFaceCentres[pointi];
1010
1011 // Bin according to surface normal
1012 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1013
1014 //- Bins of differing normals:
1015 // - one normal : flat(tish) surface
1016 // - two normals : geometric feature edge
1017 // - three normals: geometric feature point
1018 // - four normals : too complex a feature
1019 surfacePoints.clear();
1020 surfaceNormals.clear();
1021
1022 //- From face to above normals bin
1023 faceToNormalBin.setSize(pfDisp.size());
1024 faceToNormalBin = -1;
1025
1026 forAll(pfSurfNormals, i)
1027 {
1028 const point& fc = pfCentres[i];
1029 const vector& fSNormal = pfSurfNormals[i];
1030 const vector& fDisp = pfDisp[i];
1031
1032 // What to do with very far attraction? For now just ignore the face
1033 if (magSqr(fDisp) < sqr(snapDist[pointi]) && mag(fSNormal) > VSMALL)
1034 {
1035 const point pt = fc + fDisp;
1036
1037 // Do we already have surface normal?
1038 faceToNormalBin[i] = findNormal
1039 (
1040 featureCos,
1041 fSNormal,
1042 surfaceNormals
1043 );
1044
1045 if (faceToNormalBin[i] != -1)
1046 {
1047 // Same normal
1048 }
1049 else
1050 {
1051 // Now check if the planes go through the same edge or point
1052
1053 if (surfacePoints.size() <= 1)
1054 {
1055 surfacePoints.append(pt);
1056 faceToNormalBin[i] = surfaceNormals.size();
1057 surfaceNormals.append(fSNormal);
1058 }
1059 else if (surfacePoints.size() == 2)
1060 {
1061 plane pl0(surfacePoints[0], surfaceNormals[0]);
1062 plane pl1(surfacePoints[1], surfaceNormals[1]);
1063 plane::ray r(pl0.planeIntersect(pl1));
1064 vector featureNormal = r.dir() / mag(r.dir());
1065
1066 if (mag(fSNormal&featureNormal) >= 0.001)
1067 {
1068 // Definitely makes a feature point
1069 surfacePoints.append(pt);
1070 faceToNormalBin[i] = surfaceNormals.size();
1071 surfaceNormals.append(fSNormal);
1072 }
1073 }
1074 else if (surfacePoints.size() == 3)
1075 {
1076 // Have already feature point. See if this new plane is
1077 // the same point or not.
1078 plane pl0(surfacePoints[0], surfaceNormals[0]);
1079 plane pl1(surfacePoints[1], surfaceNormals[1]);
1080 plane pl2(surfacePoints[2], surfaceNormals[2]);
1081 point p012(pl0.planePlaneIntersect(pl1, pl2));
1082
1083 plane::ray r(pl0.planeIntersect(pl1));
1084 vector featureNormal = r.dir() / mag(r.dir());
1085 if (mag(fSNormal&featureNormal) >= 0.001)
1086 {
1087 plane pl3(pt, fSNormal);
1088 point p013(pl0.planePlaneIntersect(pl1, pl3));
1089
1090 if (mag(p012-p013) > snapDist[pointi])
1091 {
1092 // Different feature point
1093 surfacePoints.append(pt);
1094 faceToNormalBin[i] = surfaceNormals.size();
1095 surfaceNormals.append(fSNormal);
1096 }
1097 }
1098 }
1099 }
1100 }
1101 }
1102
1103
1104 //const point& pt = pp.localPoints()[pointi];
1105 const point& pt = ppLocalPoints[pointi];
1106
1107 // Check the number of directions
1108 if (surfaceNormals.size() == 1)
1109 {
1110 // Normal distance to plane
1111 vector d =
1112 ((surfacePoints[0]-pt) & surfaceNormals[0])
1113 *surfaceNormals[0];
1114
1115 // Trim to snap distance
1116 if (magSqr(d) > sqr(snapDist[pointi]))
1117 {
1118 d *= Foam::sqrt(sqr(snapDist[pointi])/magSqr(d));
1119 }
1120
1121 patchAttraction = d;
1122
1123 // Store constraints
1124 patchConstraint.applyConstraint(surfaceNormals[0]);
1125 }
1126 else if (surfaceNormals.size() == 2)
1127 {
1128 plane pl0(surfacePoints[0], surfaceNormals[0]);
1129 plane pl1(surfacePoints[1], surfaceNormals[1]);
1130 plane::ray r(pl0.planeIntersect(pl1));
1131 vector n = r.dir() / mag(r.dir());
1132
1133 // Get nearest point on infinite ray
1134 vector d = r.refPoint()-pt;
1135 d -= (d&n)*n;
1136
1137 // Trim to snap distance
1138 if (magSqr(d) > sqr(snapDist[pointi]))
1139 {
1140 d *= Foam::sqrt(sqr(snapDist[pointi])/magSqr(d));
1141 }
1142
1143 patchAttraction = d;
1144
1145 // Store constraints
1146 patchConstraint.applyConstraint(surfaceNormals[0]);
1147 patchConstraint.applyConstraint(surfaceNormals[1]);
1148 }
1149 else if (surfaceNormals.size() == 3)
1150 {
1151 // Calculate point from the faces.
1152 plane pl0(surfacePoints[0], surfaceNormals[0]);
1153 plane pl1(surfacePoints[1], surfaceNormals[1]);
1154 plane pl2(surfacePoints[2], surfaceNormals[2]);
1155 point cornerPt(pl0.planePlaneIntersect(pl1, pl2));
1156 vector d = cornerPt - pt;
1157
1158 // Trim to snap distance
1159 if (magSqr(d) > sqr(snapDist[pointi]))
1160 {
1161 d *= Foam::sqrt(sqr(snapDist[pointi])/magSqr(d));
1162 }
1163
1164 patchAttraction = d;
1165
1166 // Store constraints
1167 patchConstraint.applyConstraint(surfaceNormals[0]);
1168 patchConstraint.applyConstraint(surfaceNormals[1]);
1169 patchConstraint.applyConstraint(surfaceNormals[2]);
1170 }
1171}
1172
1173
1174// Special version that calculates attraction in one go
1175void Foam::snappySnapDriver::featureAttractionUsingReconstruction
1176(
1177 const label iter,
1178 const scalar featureCos,
1179
1181 const pointField& ppLocalPoints,
1182 const scalarField& snapDist,
1183 const vectorField& nearestDisp,
1184
1185 const List<List<point>>& pointFaceSurfNormals,
1186 const List<List<point>>& pointFaceDisp,
1187 const List<List<point>>& pointFaceCentres,
1188 const labelListList& pointFacePatchID,
1189
1190 vectorField& patchAttraction,
1191 List<pointConstraint>& patchConstraints
1192) const
1193{
1194 autoPtr<OBJstream> feStr;
1195 autoPtr<OBJstream> fpStr;
1197 {
1198 feStr.reset
1199 (
1200 new OBJstream
1201 (
1202 meshRefiner_.mesh().time().path()
1203 / "implicitFeatureEdge_" + name(iter) + ".obj"
1204 )
1205 );
1206 Info<< "Dumping implicit feature-edge direction to "
1207 << feStr().name() << endl;
1208
1209 fpStr.reset
1210 (
1211 new OBJstream
1212 (
1213 meshRefiner_.mesh().time().path()
1214 / "implicitFeaturePoint_" + name(iter) + ".obj"
1215 )
1216 );
1217 Info<< "Dumping implicit feature-point direction to "
1218 << fpStr().name() << endl;
1219 }
1220
1221
1222 DynamicList<point> surfacePoints(4);
1223 DynamicList<vector> surfaceNormals(4);
1224 labelList faceToNormalBin;
1225
1226 forAll(ppLocalPoints, pointi)
1227 {
1228 vector attraction = Zero;
1229 pointConstraint constraint;
1230
1231 featureAttractionUsingReconstruction
1232 (
1233 iter,
1234 featureCos,
1235
1236 pp,
1237 ppLocalPoints,
1238 snapDist,
1239 nearestDisp,
1240
1241 pointi,
1242
1243 pointFaceSurfNormals,
1244 pointFaceDisp,
1245 pointFaceCentres,
1246 pointFacePatchID,
1247
1248 surfacePoints,
1249 surfaceNormals,
1250 faceToNormalBin,
1251
1252 attraction,
1253 constraint
1254 );
1255
1256 if
1257 (
1258 (constraint.first() > patchConstraints[pointi].first())
1259 || (
1260 (constraint.first() == patchConstraints[pointi].first())
1261 && (magSqr(attraction) < magSqr(patchAttraction[pointi]))
1262 )
1263 )
1264 {
1265 patchAttraction[pointi] = attraction;
1266 patchConstraints[pointi] = constraint;
1267
1268 //const point& pt = pp.localPoints()[pointi];
1269 const point& pt = ppLocalPoints[pointi];
1270
1271 if (feStr && patchConstraints[pointi].first() == 2)
1272 {
1273 feStr().writeLine(pt, pt+patchAttraction[pointi]);
1274 }
1275 else if (fpStr && patchConstraints[pointi].first() == 3)
1276 {
1277 fpStr().writeLine(pt, pt+patchAttraction[pointi]);
1278 }
1279 }
1280 }
1281}
1282
1283
1284void Foam::snappySnapDriver::stringFeatureEdges
1285(
1286 const label iter,
1287 const scalar featureCos,
1288
1290 const pointField& ppLocalPoints,
1291 const scalarField& snapDist,
1292
1293 const vectorField& rawPatchAttraction,
1294 const List<pointConstraint>& rawPatchConstraints,
1295
1296 vectorField& patchAttraction,
1297 List<pointConstraint>& patchConstraints
1298) const
1299{
1300 // Snap edges to feature edges
1301 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1302 // Walk existing edges and snap remaining ones (that are marked as
1303 // feature edges in rawPatchConstraints)
1304
1305 // What this does is fill in any faces where not all points
1306 // on the face are being attracted:
1307 /*
1308 +
1309 / \
1310 / \
1311 ---+ +---
1312 \ /
1313 \ /
1314 +
1315 */
1316 // so the top and bottom will never get attracted since the nearest
1317 // back from the feature edge will always be one of the left or right
1318 // points since the face is diamond like. So here we walk the feature edges
1319 // and add any non-attracted points.
1320
1321
1322 while (true)
1323 {
1324 label nChanged = 0;
1325
1326 const labelListList& pointEdges = pp.pointEdges();
1327 forAll(pointEdges, pointi)
1328 {
1329 if (patchConstraints[pointi].first() == 2)
1330 {
1331 //const point& pt = pp.localPoints()[pointi];
1332 const point& pt = ppLocalPoints[pointi];
1333 const labelList& pEdges = pointEdges[pointi];
1334 const vector& featVec = patchConstraints[pointi].second();
1335
1336 // Detect whether there are edges in both directions.
1337 // (direction along the feature edge that is)
1338 bool hasPos = false;
1339 bool hasNeg = false;
1340
1341 forAll(pEdges, pEdgei)
1342 {
1343 const edge& e = pp.edges()[pEdges[pEdgei]];
1344 label nbrPointi = e.otherVertex(pointi);
1345
1346 if (patchConstraints[nbrPointi].first() > 1)
1347 {
1348 //const point& nbrPt = pp.localPoints()[nbrPointi];
1349 const point& nbrPt = ppLocalPoints[nbrPointi];
1350 const point featPt =
1351 nbrPt + patchAttraction[nbrPointi];
1352 const scalar cosAngle = (featVec & (featPt-pt));
1353
1354 if (cosAngle > 0)
1355 {
1356 hasPos = true;
1357 }
1358 else
1359 {
1360 hasNeg = true;
1361 }
1362 }
1363 }
1364
1365 if (!hasPos || !hasNeg)
1366 {
1367 //Pout<< "**Detected feature string end at "
1368 // << ppLocalPoints[pointi] << endl;
1369
1370 // No string. Assign best choice on either side
1371 label bestPosPointi = -1;
1372 scalar minPosDistSqr = GREAT;
1373 label bestNegPointi = -1;
1374 scalar minNegDistSqr = GREAT;
1375
1376 forAll(pEdges, pEdgei)
1377 {
1378 const edge& e = pp.edges()[pEdges[pEdgei]];
1379 label nbrPointi = e.otherVertex(pointi);
1380
1381 if
1382 (
1383 patchConstraints[nbrPointi].first() <= 1
1384 && rawPatchConstraints[nbrPointi].first() > 1
1385 )
1386 {
1387 const vector& nbrFeatVec =
1388 rawPatchConstraints[pointi].second();
1389
1390 if (mag(featVec&nbrFeatVec) > featureCos)
1391 {
1392 // nbrPointi attracted to sameish feature
1393 // Note: also check on position.
1394
1395 scalar d2 = magSqr
1396 (
1397 rawPatchAttraction[nbrPointi]
1398 );
1399
1400 const point featPt =
1401 //pp.localPoints()[nbrPointi]
1402 ppLocalPoints[nbrPointi]
1403 + rawPatchAttraction[nbrPointi];
1404 const scalar cosAngle =
1405 (featVec & (featPt-pt));
1406
1407 if (cosAngle > 0)
1408 {
1409 if (!hasPos && d2 < minPosDistSqr)
1410 {
1411 minPosDistSqr = d2;
1412 bestPosPointi = nbrPointi;
1413 }
1414 }
1415 else
1416 {
1417 if (!hasNeg && d2 < minNegDistSqr)
1418 {
1419 minNegDistSqr = d2;
1420 bestNegPointi = nbrPointi;
1421 }
1422 }
1423 }
1424 }
1425 }
1426
1427 if (bestPosPointi != -1)
1428 {
1429 // Use reconstructed-feature attraction. Use only
1430 // part of it since not sure...
1431 //const point& bestPt =
1432 // ppLocalPoints[bestPosPointi];
1433 //Pout<< "**Overriding point " << bestPt
1434 // << " on reconstructed feature edge at "
1435 // << rawPatchAttraction[bestPosPointi]+bestPt
1436 // << " to attracted-to-feature-edge." << endl;
1437 patchAttraction[bestPosPointi] =
1438 0.5*rawPatchAttraction[bestPosPointi];
1439 patchConstraints[bestPosPointi] =
1440 rawPatchConstraints[bestPosPointi];
1441
1442 nChanged++;
1443 }
1444 if (bestNegPointi != -1)
1445 {
1446 // Use reconstructed-feature attraction. Use only
1447 // part of it since not sure...
1448 //const point& bestPt =
1449 // ppLocalPoints[bestNegPointi];
1450 //Pout<< "**Overriding point " << bestPt
1451 // << " on reconstructed feature edge at "
1452 // << rawPatchAttraction[bestNegPointi]+bestPt
1453 // << " to attracted-to-feature-edge." << endl;
1454 patchAttraction[bestNegPointi] =
1455 0.5*rawPatchAttraction[bestNegPointi];
1456 patchConstraints[bestNegPointi] =
1457 rawPatchConstraints[bestNegPointi];
1458
1459 nChanged++;
1460 }
1461 }
1462 }
1463 }
1464
1465
1466 reduce(nChanged, sumOp<label>());
1467 Info<< "Stringing feature edges : changed " << nChanged << " points"
1468 << endl;
1469 if (nChanged == 0)
1470 {
1471 break;
1472 }
1473 }
1474}
1475
1476
1477void Foam::snappySnapDriver::releasePointsNextToMultiPatch
1478(
1479 const label iter,
1480 const scalar featureCos,
1481
1483 const pointField& ppLocalPoints,
1484 const scalarField& snapDist,
1485
1486 const List<List<point>>& pointFaceCentres,
1487 const labelListList& pointFacePatchID,
1488
1489 const vectorField& rawPatchAttraction,
1490 const List<pointConstraint>& rawPatchConstraints,
1491
1492 vectorField& patchAttraction,
1493 List<pointConstraint>& patchConstraints
1494) const
1495{
1496 autoPtr<OBJstream> multiPatchStr;
1498 {
1499 multiPatchStr.reset
1500 (
1501 new OBJstream
1502 (
1503 meshRefiner_.mesh().time().path()
1504 / "multiPatch_" + name(iter) + ".obj"
1505 )
1506 );
1507 Info<< "Dumping removed constraints due to same-face"
1508 << " multi-patch points to "
1509 << multiPatchStr().name() << endl;
1510 }
1511
1512
1513 // 1. Mark points on multiple patches
1514 bitSet isMultiPatchPoint(pp.size());
1515
1516 forAll(pointFacePatchID, pointi)
1517 {
1518 pointIndexHit multiPatchPt = findMultiPatchPoint
1519 (
1520 ppLocalPoints[pointi],
1521 pointFacePatchID[pointi],
1522 pointFaceCentres[pointi]
1523 );
1524 isMultiPatchPoint.set(pointi, multiPatchPt.hit());
1525 }
1526
1527 // 2. Make sure multi-patch points are also attracted
1528 forAll(isMultiPatchPoint, pointi)
1529 {
1530 if (isMultiPatchPoint.test(pointi))
1531 {
1532 if
1533 (
1534 patchConstraints[pointi].first() <= 1
1535 && rawPatchConstraints[pointi].first() > 1
1536 )
1537 {
1538 patchAttraction[pointi] = rawPatchAttraction[pointi];
1539 patchConstraints[pointi] = rawPatchConstraints[pointi];
1540
1541 //if (multiPatchStr)
1542 //{
1543 // Pout<< "Adding constraint on multiPatchPoint:"
1544 // << ppLocalPoints[pointi]
1545 // << " constraint:" << patchConstraints[pointi]
1546 // << " attraction:" << patchAttraction[pointi]
1547 // << endl;
1548 //}
1549 }
1550 }
1551 }
1552
1553 // Up to here it is all parallel ok.
1554
1555
1556 // 3. Knock out any attraction on faces with multi-patch points
1557 label nChanged = 0;
1558 forAll(pp.localFaces(), facei)
1559 {
1560 const face& f = pp.localFaces()[facei];
1561
1562 label nMultiPatchPoints = 0;
1563 forAll(f, fp)
1564 {
1565 label pointi = f[fp];
1566 if
1567 (
1568 isMultiPatchPoint.test(pointi)
1569 && patchConstraints[pointi].first() > 1
1570 )
1571 {
1572 ++nMultiPatchPoints;
1573 }
1574 }
1575
1576 if (nMultiPatchPoints > 0)
1577 {
1578 forAll(f, fp)
1579 {
1580 label pointi = f[fp];
1581 if
1582 (
1583 !isMultiPatchPoint.test(pointi)
1584 && patchConstraints[pointi].first() > 1
1585 )
1586 {
1587 //Pout<< "Knocking out constraint"
1588 // << " on non-multiPatchPoint:"
1589 // << ppLocalPoints[pointi] << endl;
1590 patchAttraction[pointi] = Zero;
1591 patchConstraints[pointi] = pointConstraint();
1592 nChanged++;
1593
1594 if (multiPatchStr)
1595 {
1596 multiPatchStr().write(ppLocalPoints[pointi]);
1597 }
1598 }
1599 }
1600 }
1601 }
1602
1603 reduce(nChanged, sumOp<label>());
1604 Info<< "Removing constraints near multi-patch points : changed "
1605 << nChanged << " points" << endl;
1606}
1607
1608
1609Foam::labelPair Foam::snappySnapDriver::findDiagonalAttraction
1610(
1612 const vectorField& patchAttraction,
1613 const List<pointConstraint>& patchConstraints,
1614 const label facei
1615) const
1616{
1617 const face& f = pp.localFaces()[facei];
1618 // For now just detect any attraction. Improve this to look at
1619 // actual attraction position and orientation
1620
1621 labelPair attractIndices(-1, -1);
1622
1623 if (f.size() >= 4)
1624 {
1625 for (label startFp = 0; startFp < f.size()-2; startFp++)
1626 {
1627 label minFp = f.rcIndex(startFp);
1628
1629 for
1630 (
1631 label endFp = f.fcIndex(f.fcIndex(startFp));
1632 endFp < f.size() && endFp != minFp;
1633 endFp++
1634 )
1635 {
1636 if
1637 (
1638 patchConstraints[f[startFp]].first() >= 2
1639 && patchConstraints[f[endFp]].first() >= 2
1640 )
1641 {
1642 attractIndices = labelPair(startFp, endFp);
1643 break;
1644 }
1645 }
1646 }
1647 }
1648 return attractIndices;
1649}
1650
1651
1652bool Foam::snappySnapDriver::isSplitAlignedWithFeature
1653(
1654 const scalar featureCos,
1655 const point& p0,
1656 const pointConstraint& pc0,
1657 const point& p1,
1658 const pointConstraint& pc1
1659) const
1660{
1661 vector d(p1-p0);
1662 scalar magD = mag(d);
1663 if (magD < VSMALL)
1664 {
1665 // Two diagonal points already colocated?
1666 return false;
1667 }
1668 else
1669 {
1670 d /= magD;
1671
1672 // Is diagonal d aligned with at least one of the feature
1673 // edges?
1674
1675 if (pc0.first() == 2 && mag(d & pc0.second()) > featureCos)
1676 {
1677 return true;
1678 }
1679 else if (pc1.first() == 2 && mag(d & pc1.second()) > featureCos)
1680 {
1681 return true;
1682 }
1683 else
1684 {
1685 return false;
1686 }
1687 }
1688}
1689
1690
1691// Is situation very concave
1692bool Foam::snappySnapDriver::isConcave
1693(
1694 const point& c0,
1695 const vector& area0,
1696 const point& c1,
1697 const vector& area1,
1698 const scalar concaveCos
1699) const
1700{
1701 vector n0 = area0;
1702 scalar magN0 = mag(n0);
1703 if (magN0 < VSMALL)
1704 {
1705 // Zero area face. What to return? For now disable splitting.
1706 return true;
1707 }
1708 n0 /= magN0;
1709
1710 // Distance from c1 to plane of face0
1711 scalar d = (c1-c0)&n0;
1712
1713 if (d <= 0)
1714 {
1715 // Convex (face1 centre on 'inside' of face0)
1716 return false;
1717 }
1718 else
1719 {
1720 // Is a bit or very concave?
1721 vector n1 = area1;
1722 scalar magN1 = mag(n1);
1723 if (magN1 < VSMALL)
1724 {
1725 // Zero area face. See above
1726 return true;
1727 }
1728 n1 /= magN1;
1729
1730 if ((n0&n1) < concaveCos)
1731 {
1732 return true;
1733 }
1734 else
1735 {
1736 return false;
1737 }
1738 }
1739}
1740
1741
1742Foam::labelPair Foam::snappySnapDriver::findDiagonalAttraction
1743(
1744 const scalar featureCos,
1745 const scalar concaveCos,
1746 const scalar minAreaRatio,
1748 const pointField& ppLocalPoints,
1749 const vectorField& patchAttr,
1750 const List<pointConstraint>& patchConstraints,
1751 const vectorField& nearestAttr,
1752 const vectorField& nearestNormal,
1753 const label facei,
1754
1756 DynamicField<point>& points1
1757) const
1758{
1759 const face& localF = pp.localFaces()[facei];
1760
1761 labelPair attractIndices(-1, -1);
1762
1763 if (localF.size() >= 4)
1764 {
1768 //const polyMesh& mesh = meshRefiner_.mesh();
1769 //label meshFacei = pp.addressing()[facei];
1770 //const face& meshF = mesh.faces()[meshFacei];
1771 //label celli = mesh.faceOwner()[meshFacei];
1772 //const labelList& cPoints = mesh.cellPoints(celli);
1773 //
1774 //point cc(mesh.points()[meshF[0]]);
1775 //for (label i = 1; i < meshF.size(); i++)
1776 //{
1777 // cc += mesh.points()[meshF[i]]+patchAttr[localF[i]];
1778 //}
1779 //forAll(cPoints, i)
1780 //{
1781 // label pointi = cPoints[i];
1782 // if (!meshF.found(pointi))
1783 // {
1784 // cc += mesh.points()[pointi];
1785 // }
1786 //}
1787 //cc /= cPoints.size();
1789 //
1790 //const scalar vol = pyrVol(pp, patchAttr, localF, cc);
1791 //const scalar area = localF.mag(localPts);
1792
1793
1794
1795 // Try all diagonal cuts
1796 // ~~~~~~~~~~~~~~~~~~~~~
1797
1798 face f0(3);
1799 face f1(3);
1800
1801 for (label startFp = 0; startFp < localF.size()-2; startFp++)
1802 {
1803 label minFp = localF.rcIndex(startFp);
1804
1805 for
1806 (
1807 label endFp = localF.fcIndex(localF.fcIndex(startFp));
1808 endFp < localF.size() && endFp != minFp;
1809 endFp++
1810 )
1811 {
1812 label startPti = localF[startFp];
1813 label endPti = localF[endFp];
1814
1815 const pointConstraint& startPc = patchConstraints[startPti];
1816 const pointConstraint& endPc = patchConstraints[endPti];
1817
1818 if (startPc.first() >= 2 && endPc.first() >= 2)
1819 {
1820 if (startPc.first() == 2 || endPc.first() == 2)
1821 {
1822 // Check if
1823 // - sameish feature edge normal
1824 // - diagonal aligned with feature edge normal
1825 point start =
1826 ppLocalPoints[startPti]+patchAttr[startPti];
1827 point end =
1828 ppLocalPoints[endPti]+patchAttr[endPti];
1829
1830 if
1831 (
1832 !isSplitAlignedWithFeature
1833 (
1834 featureCos,
1835 start,
1836 startPc,
1837 end,
1838 endPc
1839 )
1840 )
1841 {
1842 // Attract to different features. No need to
1843 // introduce split
1844 continue;
1845 }
1846 }
1847
1848
1849
1850 // Form two faces
1851 // ~~~~~~~~~~~~~~
1852 // Predict position of faces. End points of the faces
1853 // attract to the feature
1854 // and all the other points just attract to the nearest
1855
1856 // face0
1857
1858 f0.setSize(endFp-startFp+1);
1859 label i0 = 0;
1860 for (label fp = startFp; fp <= endFp; fp++)
1861 {
1862 f0[i0++] = localF[fp];
1863 }
1864
1865 // Get compact face and points
1866 const face compact0(identity(f0.size()));
1867 points0.clear();
1868 points0.append(ppLocalPoints[f0[0]] + patchAttr[f0[0]]);
1869 for (label fp=1; fp < f0.size()-1; fp++)
1870 {
1871 label pi = f0[fp];
1872 points0.append(ppLocalPoints[pi] + nearestAttr[pi]);
1873 }
1875 (
1876 ppLocalPoints[f0.last()] + patchAttr[f0.last()]
1877 );
1878
1879
1880 // face1
1881
1882 f1.setSize(localF.size()+2-f0.size());
1883 label i1 = 0;
1884 for
1885 (
1886 label fp = endFp;
1887 fp != startFp;
1888 fp = localF.fcIndex(fp)
1889 )
1890 {
1891 f1[i1++] = localF[fp];
1892 }
1893 f1[i1++] = localF[startFp];
1894
1895
1896 // Get compact face and points
1897 const face compact1(identity(f1.size()));
1898 points1.clear();
1899 points1.append(ppLocalPoints[f1[0]] + patchAttr[f1[0]]);
1900 for (label fp=1; fp < f1.size()-1; fp++)
1901 {
1902 label pi = f1[fp];
1903 points1.append(ppLocalPoints[pi] + nearestAttr[pi]);
1904 }
1905 points1.append
1906 (
1907 ppLocalPoints[f1.last()] + patchAttr[f1.last()]
1908 );
1909
1910
1911
1912 // Avoid splitting concave faces
1913 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1914
1915 if
1916 (
1917 isConcave
1918 (
1919 compact0.centre(points0),
1920 compact0.areaNormal(points0),
1921 compact1.centre(points1),
1922 compact1.areaNormal(points1),
1923 concaveCos
1924 )
1925 )
1926 {
1928 //Pout<< "# f0:" << f0 << endl;
1929 //forAll(p, i)
1930 //{
1931 // meshTools::writeOBJ(Pout, points0[i]);
1932 //}
1933 //Pout<< "# f1:" << f1 << endl;
1934 //forAll(p, i)
1935 //{
1936 // meshTools::writeOBJ(Pout, points1[i]);
1937 //}
1938 }
1939 else
1940 {
1941 // Existing areas
1942 const scalar area0 = f0.mag(ppLocalPoints);
1943 const scalar area1 = f1.mag(ppLocalPoints);
1944
1945 if
1946 (
1947 area0/area1 >= minAreaRatio
1948 && area1/area0 >= minAreaRatio
1949 )
1950 {
1951 attractIndices = labelPair(startFp, endFp);
1952 }
1953 }
1954 }
1955 }
1956 }
1957 }
1958 return attractIndices;
1959}
1960
1961
1962void Foam::snappySnapDriver::splitDiagonals
1963(
1964 const scalar featureCos,
1965 const scalar concaveCos,
1966 const scalar minAreaRatio,
1967
1969 const pointField& ppLocalPoints,
1970 const vectorField& nearestAttraction,
1971 const vectorField& nearestNormal,
1972 const List<labelList>& pointFacePatchID,
1973
1974 vectorField& patchAttraction,
1975 List<pointConstraint>& patchConstraints,
1976
1977 DynamicList<label>& splitFaces,
1978 DynamicList<labelPair>& splits,
1979 DynamicList<labelPair>& splitPatches
1980) const
1981{
1982 const labelList& bFaces = pp.addressing();
1983
1984 splitFaces.clear();
1985 splitFaces.setCapacity(bFaces.size());
1986 splits.clear();
1987 splits.setCapacity(bFaces.size());
1988 splitPatches.clear();
1989 splitPatches.setCapacity(bFaces.size());
1990
1991 // Work arrays for storing points of face
1992 DynamicField<point> facePoints0;
1993 DynamicField<point> facePoints1;
1994 face f0, f1;
1995
1996 forAll(bFaces, facei)
1997 {
1998 const labelPair split
1999 (
2000 findDiagonalAttraction
2001 (
2002 featureCos,
2003 concaveCos,
2004 minAreaRatio,
2005
2006 pp,
2007 ppLocalPoints,
2008 patchAttraction,
2009 patchConstraints,
2010
2011 nearestAttraction,
2012 nearestNormal,
2013 facei,
2014
2015 facePoints0,
2016 facePoints1
2017 )
2018 );
2019
2020 if (split != labelPair(-1, -1))
2021 {
2022 splitFaces.append(bFaces[facei]);
2023 splits.append(split);
2024 splitPatches.append(labelPair(-1, -1));
2025
2026 const face& f = pp.localFaces()[facei];
2027
2028 // Knock out other attractions on face
2029 forAll(f, fp)
2030 {
2031 if
2032 (
2033 fp != split[0]
2034 && fp != split[1]
2035 && patchConstraints[f[fp]].first() >= 2
2036 )
2037 {
2038 patchConstraints[f[fp]] = pointConstraint
2039 (
2041 (
2042 1,
2043 nearestNormal[f[fp]]
2044 )
2045 );
2046 patchAttraction[f[fp]] = nearestAttraction[f[fp]];
2047 }
2048 }
2049
2050
2051 // Detect any patches to give to the two faces
2052 labelPair& twoPatches = splitPatches.last();
2053
2054 // Split into two faces:
2055 // - face0 is from split[0] up to split[1]
2056 // - face1 is from split[1] up to split[0]
2058
2059 // Bit heuristic:
2060 // - find any inbetween vertex (i.e. not on diagonal) that
2061 // has single patch only
2062 // - if we find such a point on both sides we use these patches
2063 // - if only found a point on one side check the other side
2064 // for any other patch
2065
2066 // Find single patch
2067 for (label fp = 1; fp < f0.size()-1; fp++)
2068 {
2069 const auto& patches = pointFacePatchID[f0[fp]];
2070 if (patches.uniform())
2071 {
2072 twoPatches[0] = patches[0];
2073 }
2074 }
2075 for (label fp = 1; fp < f1.size()-1; fp++)
2076 {
2077 const auto& patches = pointFacePatchID[f1[fp]];
2078 if (patches.uniform())
2079 {
2080 twoPatches[1] = patches[0];
2081 }
2082 }
2083
2084 if (twoPatches[0] == -1)
2085 {
2086 // Find any patch on any point that differs from twoPatches[1]
2087 // (can be -1 or a unform patch)
2088
2089 for (label fp = 1; fp < f0.size()-1; fp++)
2090 {
2091 const auto& patches = pointFacePatchID[f0[fp]];
2092
2093 // Pick any patch which is not the unique patch of the
2094 // other face. TBD: count occurences and choose most used?
2095 for (const label patchi : patches)
2096 {
2097 if (patchi != twoPatches[1])
2098 {
2099 twoPatches[0] = patchi;
2100 break;
2101 }
2102 }
2103 }
2104 }
2105
2106 if (twoPatches[1] == -1)
2107 {
2108 // Find any patch on any point that differs from twoPatches[0]
2109 // (can be -1 or a unform patch)
2110
2111 for (label fp = 1; fp < f1.size()-1; fp++)
2112 {
2113 const auto& patches = pointFacePatchID[f1[fp]];
2114
2115 // Pick any patch which is not the unique patch of the
2116 // other face. TBD: count occurences and choose most used?
2117 for (const label patchi : patches)
2118 {
2119 if (patchi != twoPatches[0])
2120 {
2121 twoPatches[1] = patchi;
2122 break;
2123 }
2124 }
2125 }
2126 }
2127
2128 // Fall-back : choose any connected patch
2129 if (twoPatches[0] == -1)
2130 {
2131 twoPatches[0] = pointFacePatchID[f0[1]][0];
2132 }
2133 if (twoPatches[1] == -1)
2134 {
2135 twoPatches[1] = pointFacePatchID[f1[1]][0];
2136 }
2137
2138 //if (twoPatches != labelPair(-1, -1))
2139 //{
2140 // Pout<< "Face:" << bFaces[facei]
2141 // << " at:" << pp.faceCentres()[facei]
2142 // << " split into f0:" << flatOutput(f0)
2143 // << " patch:" << twoPatches[0]
2144 // << " and f1:" << flatOutput(f1)
2145 // << " patch:" << twoPatches[1]
2146 // << endl;
2147 //}
2148 }
2149 }
2150}
2151
2152
2153void Foam::snappySnapDriver::avoidDiagonalAttraction
2154(
2155 const label iter,
2156 const scalar featureCos,
2157
2159 const pointField& ppLocalPoints,
2160
2161 vectorField& patchAttraction,
2162 List<pointConstraint>& patchConstraints
2163) const
2164{
2165 forAll(pp.localFaces(), facei)
2166 {
2167 const face& f = pp.localFaces()[facei];
2168
2169 labelPair diag = findDiagonalAttraction
2170 (
2171 pp,
2172 patchAttraction,
2173 patchConstraints,
2174 facei
2175 );
2176
2177 if (diag[0] != -1 && diag[1] != -1)
2178 {
2179 // Found two diagonal points that being attracted.
2180 // For now just attract my one to the average of those.
2181 const label i0 = f[diag[0]];
2182 const point pt0 =
2183 ppLocalPoints[i0]+patchAttraction[i0];
2184 const label i1 = f[diag[1]];
2185 const point pt1 =
2186 ppLocalPoints[i1]+patchAttraction[i1];
2187 const point mid = 0.5*(pt0+pt1);
2188
2189 const scalar cosAngle = mag
2190 (
2191 patchConstraints[i0].second()
2192 & patchConstraints[i1].second()
2193 );
2194
2195 //Pout<< "Found diagonal attraction at indices:"
2196 // << diag[0]
2197 // << " and " << diag[1]
2198 // << " with cosAngle:" << cosAngle
2199 // << " mid:" << mid << endl;
2200
2201 if (cosAngle > featureCos)
2202 {
2203 // The two feature edges are roughly in the same direction.
2204 // Add the nearest of the other points in the face as
2205 // attractor
2206 label minFp = -1;
2207 scalar minDistSqr = GREAT;
2208 forAll(f, fp)
2209 {
2210 label pointi = f[fp];
2211 if (patchConstraints[pointi].first() <= 1)
2212 {
2213 scalar distSqr = mid.distSqr(ppLocalPoints[pointi]);
2214 if (distSqr < minDistSqr)
2215 {
2216 minFp = fp;
2217 }
2218 }
2219 }
2220 if (minFp != -1)
2221 {
2222 label minPointi = f[minFp];
2223 patchAttraction[minPointi] =
2224 mid-ppLocalPoints[minPointi];
2225 patchConstraints[minPointi] = patchConstraints[f[diag[0]]];
2226 }
2227 }
2228 else
2229 {
2230 //Pout<< "Diagonal attractors at" << nl
2231 // << " pt0:" << pt0
2232 // << " constraint:"
2233 // << patchConstraints[i0].second() << nl
2234 // << " pt1:" << pt1
2235 // << " constraint:"
2236 // << patchConstraints[i1].second() << nl
2237 // << " make too large an angle:"
2238 // << mag
2239 // (
2240 // patchConstraints[i0].second()
2241 // & patchConstraints[i1].second()
2242 // )
2243 // << endl;
2244 }
2245 }
2246 }
2247}
2248
2249
2250Foam::Tuple2<Foam::label, Foam::pointIndexHit>
2251Foam::snappySnapDriver::findNearFeatureEdge
2252(
2253 const bool isRegionEdge,
2254
2256 const pointField& ppLocalPoints,
2257 const scalarField& snapDist,
2258 const label pointi,
2259 const point& estimatedPt,
2260
2261 List<List<DynamicList<point>>>& edgeAttractors,
2262 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
2263 vectorField& patchAttraction,
2264 List<pointConstraint>& patchConstraints
2265) const
2266{
2267 const refinementFeatures& features = meshRefiner_.features();
2268
2269 labelList nearEdgeFeat;
2270 List<pointIndexHit> nearEdgeInfo;
2271 vectorField nearNormal;
2272
2273 if (isRegionEdge)
2274 {
2275 features.findNearestRegionEdge
2276 (
2277 pointField(1, estimatedPt),
2278 scalarField(1, sqr(snapDist[pointi])),
2279 nearEdgeFeat,
2280 nearEdgeInfo,
2281 nearNormal
2282 );
2283 }
2284 else
2285 {
2286 features.findNearestEdge
2287 (
2288 pointField(1, estimatedPt),
2289 scalarField(1, sqr(snapDist[pointi])),
2290 nearEdgeFeat,
2291 nearEdgeInfo,
2292 nearNormal
2293 );
2294 }
2295
2296 const pointIndexHit& nearInfo = nearEdgeInfo[0];
2297 label feati = nearEdgeFeat[0];
2298
2299 if (nearInfo.hit())
2300 {
2301 // So we have a point on the feature edge. Use this
2302 // instead of our estimate from planes.
2303 edgeAttractors[feati][nearInfo.index()].append(nearInfo.point());
2304 pointConstraint c(Tuple2<label, vector>(2, nearNormal[0]));
2305 edgeConstraints[feati][nearInfo.index()].append(c);
2306
2307 // Store for later use
2308 patchAttraction[pointi] = nearInfo.point()-ppLocalPoints[pointi];
2309 patchConstraints[pointi] = c;
2310 }
2311 return Tuple2<label, pointIndexHit>(feati, nearInfo);
2312}
2313
2314
2315Foam::Tuple2<Foam::label, Foam::pointIndexHit>
2316Foam::snappySnapDriver::findNearFeaturePoint
2317(
2318 const bool isRegionPoint,
2319
2321 const pointField& ppLocalPoints,
2322 const scalarField& snapDist,
2323 const label pointi,
2324 const point& estimatedPt,
2325
2326 // Feature-point to pp point
2327 List<labelList>& pointAttractor,
2329 // Feature-edge to pp point
2330 List<List<DynamicList<point>>>& edgeAttractors,
2331 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
2332 // pp point to nearest feature
2333 vectorField& patchAttraction,
2334 List<pointConstraint>& patchConstraints
2335) const
2336{
2337 const refinementFeatures& features = meshRefiner_.features();
2338
2339 // Search for for featurePoints only! This ignores any non-feature points.
2340
2341 labelList nearFeat;
2342 List<pointIndexHit> nearInfo;
2343 features.findNearestPoint
2344 (
2345 pointField(1, estimatedPt),
2346 scalarField(1, sqr(snapDist[pointi])),
2347 nearFeat,
2348 nearInfo
2349 );
2350
2351 label feati = nearFeat[0];
2352
2353 if (feati != -1)
2354 {
2355 const point& pt = ppLocalPoints[pointi];
2356
2357 label featPointi = nearInfo[0].index();
2358 const point& featPt = nearInfo[0].hitPoint();
2359 scalar distSqr = featPt.distSqr(pt);
2360
2361 // Check if already attracted
2362 label oldPointi = pointAttractor[feati][featPointi];
2363
2364 if (oldPointi != -1)
2365 {
2366 // Check distance
2367 if (distSqr >= featPt.distSqr(ppLocalPoints[oldPointi]))
2368 {
2369 // oldPointi nearest. Keep.
2370 feati = -1;
2371 featPointi = -1;
2372 }
2373 else
2374 {
2375 // Current pointi nearer.
2376 pointAttractor[feati][featPointi] = pointi;
2377 pointConstraints[feati][featPointi].first() = 3;
2378 pointConstraints[feati][featPointi].second() = Zero;
2379
2380 // Store for later use
2381 patchAttraction[pointi] = featPt-pt;
2382 patchConstraints[pointi] = pointConstraints[feati][featPointi];
2383
2384 // Reset oldPointi to nearest on feature edge
2385 patchAttraction[oldPointi] = Zero;
2386 patchConstraints[oldPointi] = pointConstraint();
2387
2388 findNearFeatureEdge
2389 (
2390 isRegionPoint, // search region edges only
2391
2392 pp,
2393 ppLocalPoints,
2394 snapDist,
2395 oldPointi,
2396 ppLocalPoints[oldPointi],
2397
2398 edgeAttractors,
2399 edgeConstraints,
2400 patchAttraction,
2401 patchConstraints
2402 );
2403 }
2404 }
2405 else
2406 {
2407 // Current pointi nearer.
2408 pointAttractor[feati][featPointi] = pointi;
2409 pointConstraints[feati][featPointi].first() = 3;
2410 pointConstraints[feati][featPointi].second() = Zero;
2411
2412 // Store for later use
2413 patchAttraction[pointi] = featPt-pt;
2414 patchConstraints[pointi] = pointConstraints[feati][featPointi];
2415 }
2416 }
2417
2418 return Tuple2<label, pointIndexHit>(feati, nearInfo[0]);
2419}
2420
2421
2422// Determines for every pp point - that is on multiple faces that form
2423// a feature - the nearest feature edge/point.
2424void Foam::snappySnapDriver::determineFeatures
2425(
2426 const label iter,
2427 const scalar featureCos,
2428 const bool multiRegionFeatureSnap,
2429 const bool strictRegionFeatureSnap, // special >=3 patch points handling
2430
2432 const pointField& ppLocalPoints,
2433 const scalarField& snapDist,
2434 const vectorField& nearestDisp,
2435
2436 const List<List<point>>& pointFaceSurfNormals,
2437 const List<List<point>>& pointFaceDisp,
2438 const List<List<point>>& pointFaceCentres,
2439 const labelListList& pointFacePatchID,
2440
2441 // Feature-point to pp point
2442 List<labelList>& pointAttractor,
2444 // Feature-edge to pp point
2445 List<List<DynamicList<point>>>& edgeAttractors,
2446 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
2447 // pp point to nearest feature
2448 vectorField& patchAttraction,
2449 List<pointConstraint>& patchConstraints
2450) const
2451{
2452 autoPtr<OBJstream> featureEdgeStr;
2453 autoPtr<OBJstream> missedEdgeStr;
2454 autoPtr<OBJstream> featurePointStr;
2455 autoPtr<OBJstream> missedMP0Str;
2456 autoPtr<OBJstream> missedMP1Str;
2457
2459 {
2460 featureEdgeStr.reset
2461 (
2462 new OBJstream
2463 (
2464 meshRefiner_.mesh().time().path()
2465 / "featureEdge_" + name(iter) + ".obj"
2466 )
2467 );
2468 Info<< "Dumping feature-edge sampling to "
2469 << featureEdgeStr().name() << endl;
2470
2471 missedEdgeStr.reset
2472 (
2473 new OBJstream
2474 (
2475 meshRefiner_.mesh().time().path()
2476 / "missedFeatureEdge_" + name(iter) + ".obj"
2477 )
2478 );
2479 Info<< "Dumping feature-edges that are too far away to "
2480 << missedEdgeStr().name() << endl;
2481
2482 featurePointStr.reset
2483 (
2484 new OBJstream
2485 (
2486 meshRefiner_.mesh().time().path()
2487 / "featurePoint_" + name(iter) + ".obj"
2488 )
2489 );
2490 Info<< "Dumping feature-point sampling to "
2491 << featurePointStr().name() << endl;
2492
2493 missedMP0Str.reset
2494 (
2495 new OBJstream
2496 (
2497 meshRefiner_.mesh().time().path()
2498 / "missedFeatureEdgeFromMPEdge_" + name(iter) + ".obj"
2499 )
2500 );
2501 Info<< "Dumping region-edges that are too far away to "
2502 << missedMP0Str().name() << endl;
2503
2504 missedMP1Str.reset
2505 (
2506 new OBJstream
2507 (
2508 meshRefiner_.mesh().time().path()
2509 / "missedFeatureEdgeFromMPPoint_" + name(iter) + ".obj"
2510 )
2511 );
2512 Info<< "Dumping region-points that are too far away to "
2513 << missedMP1Str().name() << endl;
2514 }
2515
2516
2517 DynamicList<point> surfacePoints(4);
2518 DynamicList<vector> surfaceNormals(4);
2519 labelList faceToNormalBin;
2520
2521 forAll(ppLocalPoints, pointi)
2522 {
2523 const point& pt = ppLocalPoints[pointi];
2524
2525
2526 // Determine the geometric planes the point is (approximately) on.
2527 // This is returned as a
2528 // - attraction vector
2529 // - and a constraint
2530 // (1: attract to surface, constraint is normal of plane
2531 // 2: attract to feature line, constraint is feature line direction
2532 // 3: attract to feature point, constraint is zero)
2533
2534 vector attraction = Zero;
2535 pointConstraint constraint;
2536
2537 featureAttractionUsingReconstruction
2538 (
2539 iter,
2540 featureCos,
2541
2542 pp,
2543 ppLocalPoints,
2544 snapDist,
2545 nearestDisp,
2546
2547 pointi,
2548
2549 pointFaceSurfNormals,
2550 pointFaceDisp,
2551 pointFaceCentres,
2552 pointFacePatchID, // currently not used
2553
2554 surfacePoints,
2555 surfaceNormals,
2556 faceToNormalBin,
2557
2558 attraction,
2559 constraint
2560 );
2561
2562
2563 if (strictRegionFeatureSnap)
2564 {
2565 // Bit tricky: if a point is on more than 2 points and only
2566 // attracted-to-surface upgrade it to attract-to-edge so it follows
2567 // more complex logic below. This is easier than complicating
2568 // the nearestFeatureEdge/Point below.
2569
2570 if (constraint.first() == 1)
2571 {
2572 const auto& patches = pointFacePatchID[pointi];
2573 label patch1 = -1;
2574 for (label i = 1; i < patches.size(); i++)
2575 {
2576 if (patches[i] != patches[0])
2577 {
2578 if (patch1 == -1)
2579 {
2580 patch1 = patches[i];
2581 }
2582 else if (patches[i] != patch1)
2583 {
2584 // >= on 2 patches. Upgrade to edge-attract
2585 constraint.first() = 2;
2586 //Pout<< "** 3 patch point:" << pointi
2587 // << " at:" << pt
2588 // << " patches:" << flatOutput(patches)
2589 // << " upgraded to feature-edge"
2590 // << " attracted to " << attraction << endl;
2591 break;
2592 }
2593 }
2594 }
2595 }
2596 }
2597
2598
2599 // Now combine the reconstruction with the current state of the
2600 // point. The logic is quite complicated:
2601 // - the new constraint (from reconstruction) will only win if
2602 // - the constraint is higher (feature-point wins from feature-edge
2603 // etc.)
2604 // - or the constraint is the same but the attraction distance is less
2605 //
2606 // - then this will be combined with explicit searching on the
2607 // features and optionally the analysis of the patches using the
2608 // point. This analysis can do three thing:
2609 // - the point is not on multiple patches
2610 // - the point is on multiple patches but these are also
2611 // different planes (so the region feature is also a geometric
2612 // feature)
2613 // - the point is on multiple patches some of which are on
2614 // the same plane. This is the problem one - do we assume it is
2615 // an additional constraint (feat edge upgraded to region point,
2616 // see below)?
2617 //
2618 // Reconstruction MultiRegionFeatureSnap Attraction
2619 // ------- ---------------------- -----------
2620 // surface false surface
2621 // surface true region edge
2622 // feat edge false feat edge
2623 // feat edge true and no planar regions feat edge
2624 // feat edge true and yes planar regions region point
2625 // feat point false feat point
2626 // feat point true region point
2627
2628
2629 if
2630 (
2631 (constraint.first() > patchConstraints[pointi].first())
2632 || (
2633 (constraint.first() == patchConstraints[pointi].first())
2634 && (magSqr(attraction) < magSqr(patchAttraction[pointi]))
2635 )
2636 )
2637 {
2638 patchAttraction[pointi] = attraction;
2639 patchConstraints[pointi] = constraint;
2640
2641 // Check the number of directions
2642 if (patchConstraints[pointi].first() == 1)
2643 {
2644 // Flat surface. Check for different patchIDs
2645 if (multiRegionFeatureSnap)
2646 {
2647 const point estimatedPt(pt + nearestDisp[pointi]);
2648 pointIndexHit multiPatchPt
2649 (
2650 findMultiPatchPoint
2651 (
2652 estimatedPt,
2653 pointFacePatchID[pointi],
2654 surfaceNormals,
2655 faceToNormalBin
2656 )
2657 );
2658
2659 if (multiPatchPt.hit())
2660 {
2661 // Behave like when having two surface normals so
2662 // attract to nearest feature edge (with a guess for
2663 // the multipatch point as starting point)
2665 findNearFeatureEdge
2666 (
2667 true, // isRegionEdge
2668 pp,
2669 ppLocalPoints,
2670 snapDist,
2671 pointi,
2672 multiPatchPt.point(), // estimatedPt
2673
2674 edgeAttractors,
2675 edgeConstraints,
2676
2677 patchAttraction,
2678 patchConstraints
2679 );
2680
2681 const pointIndexHit& info = nearInfo.second();
2682 if (info.hit())
2683 {
2684 // Dump
2685 if (featureEdgeStr)
2686 {
2687 featureEdgeStr().writeLine
2688 (
2689 pt,
2690 info.point()
2691 );
2692 }
2693 }
2694 else
2695 {
2696 if (missedEdgeStr)
2697 {
2698 missedEdgeStr().writeLine
2699 (
2700 pt,
2701 multiPatchPt.point()
2702 );
2703 }
2704 }
2705 }
2706 }
2707 }
2708 else if (patchConstraints[pointi].first() == 2)
2709 {
2710 // Mark point on the nearest feature edge. Note that we
2711 // only search within the surrounding since the plane
2712 // reconstruction might find a feature where there isn't one.
2713 const point estimatedPt(pt + patchAttraction[pointi]);
2714
2716
2717 // Geometric feature edge. Check for different patchIDs
2718 bool hasSnapped = false;
2719 if (multiRegionFeatureSnap)
2720 {
2721 pointIndexHit multiPatchPt
2722 (
2723 findMultiPatchPoint
2724 (
2725 estimatedPt,
2726 pointFacePatchID[pointi],
2727 surfaceNormals,
2728 faceToNormalBin
2729 )
2730 );
2731 if (multiPatchPt.hit())
2732 {
2733 if (multiPatchPt.index() == 0)
2734 {
2735 // Region edge is also a geometric feature edge
2736 nearInfo = findNearFeatureEdge
2737 (
2738 true, // isRegionEdge
2739 pp,
2740 ppLocalPoints,
2741 snapDist,
2742 pointi,
2743 estimatedPt,
2744
2745 edgeAttractors,
2746 edgeConstraints,
2747
2748 patchAttraction,
2749 patchConstraints
2750 );
2751 hasSnapped = true;
2752
2753 // Debug: dump missed feature point
2754 if (missedMP0Str && !nearInfo.second().hit())
2755 {
2756 missedMP0Str().writeLine(pt, estimatedPt);
2757 }
2758 }
2759 else
2760 {
2761 // One of planes of feature contains multiple
2762 // regions. We assume (contentious!) that the
2763 // separation between
2764 // the regions is not aligned with the geometric
2765 // feature so is an additional constraint on the
2766 // point -> is region-feature-point.
2767 nearInfo = findNearFeaturePoint
2768 (
2769 true, // isRegionPoint
2770 pp,
2771 ppLocalPoints,
2772 snapDist,
2773 pointi,
2774 estimatedPt,
2775
2776 // Feature-point to pp point
2777 pointAttractor,
2779 // Feature-edge to pp point
2780 edgeAttractors,
2781 edgeConstraints,
2782 // pp point to nearest feature
2783 patchAttraction,
2784 patchConstraints
2785 );
2786 hasSnapped = true;
2787
2788 // More contentious: if we don't find
2789 // a near feature point we will never find the
2790 // attraction to a feature edge either since
2791 // the edgeAttractors/edgeConstraints do not get
2792 // filled and we're using reverse attraction
2793 // Note that we're in multiRegionFeatureSnap which
2794 // where findMultiPatchPoint can decide the
2795 // wrong thing. So: if failed finding a near
2796 // feature point try for a feature edge
2797 if (!nearInfo.second().hit())
2798 {
2799 nearInfo = findNearFeatureEdge
2800 (
2801 true, // isRegionEdge
2802 pp,
2803 ppLocalPoints,
2804 snapDist,
2805 pointi,
2806 estimatedPt,
2807
2808 // Feature-edge to pp point
2809 edgeAttractors,
2810 edgeConstraints,
2811 // pp point to nearest feature
2812 patchAttraction,
2813 patchConstraints
2814 );
2815 }
2816
2817 // Debug: dump missed feature point
2818 if (missedMP1Str && !nearInfo.second().hit())
2819 {
2820 missedMP1Str().writeLine(pt, estimatedPt);
2821 }
2822 }
2823 }
2824 }
2825
2826 if (!hasSnapped)
2827 {
2828 // Determine nearest point on feature edge. Store
2829 // constraint
2830 // (calculated from feature edge, alternative would be to
2831 // use constraint calculated from both surfaceNormals)
2832 nearInfo = findNearFeatureEdge
2833 (
2834 false, // isRegionPoint
2835 pp,
2836 ppLocalPoints,
2837 snapDist,
2838 pointi,
2839 estimatedPt,
2840
2841 edgeAttractors,
2842 edgeConstraints,
2843
2844 patchAttraction,
2845 patchConstraints
2846 );
2847 hasSnapped = true;
2848 }
2849
2850 // Dump to obj
2851 const pointIndexHit& info = nearInfo.second();
2852 if (info.hit())
2853 {
2854 if
2855 (
2856 featurePointStr
2857 && patchConstraints[pointi].first() == 3
2858 )
2859 {
2860 featurePointStr().writeLine(pt, info.point());
2861 }
2862 else if
2863 (
2864 featureEdgeStr
2865 && patchConstraints[pointi].first() == 2
2866 )
2867 {
2868 featureEdgeStr().writeLine(pt, info.point());
2869 }
2870 }
2871 else
2872 {
2873 if (missedEdgeStr)
2874 {
2875 missedEdgeStr().writeLine(pt, estimatedPt);
2876 }
2877 }
2878 }
2879 else if (patchConstraints[pointi].first() == 3)
2880 {
2881 // Mark point on the nearest feature point.
2882 const point estimatedPt(pt + patchAttraction[pointi]);
2883
2885
2886 if (multiRegionFeatureSnap)
2887 {
2888 pointIndexHit multiPatchPt
2889 (
2890 findMultiPatchPoint
2891 (
2892 estimatedPt,
2893 pointFacePatchID[pointi],
2894 surfaceNormals,
2895 faceToNormalBin
2896 )
2897 );
2898 if (multiPatchPt.hit())
2899 {
2900 // Multiple regions
2901 nearInfo = findNearFeaturePoint
2902 (
2903 true, // isRegionPoint
2904 pp,
2905 ppLocalPoints,
2906 snapDist,
2907 pointi,
2908 estimatedPt,
2909
2910 // Feature-point to pp point
2911 pointAttractor,
2913 // Feature-edge to pp point
2914 edgeAttractors,
2915 edgeConstraints,
2916 // pp point to nearest feature
2917 patchAttraction,
2918 patchConstraints
2919 );
2920 }
2921 else
2922 {
2923 nearInfo = findNearFeaturePoint
2924 (
2925 false, // isRegionPoint
2926 pp,
2927 ppLocalPoints,
2928 snapDist,
2929 pointi,
2930 estimatedPt,
2931
2932 // Feature-point to pp point
2933 pointAttractor,
2935 // Feature-edge to pp point
2936 edgeAttractors,
2937 edgeConstraints,
2938 // pp point to nearest feature
2939 patchAttraction,
2940 patchConstraints
2941 );
2942 }
2943 }
2944 else
2945 {
2946 // No multi-patch snapping
2947 nearInfo = findNearFeaturePoint
2948 (
2949 false, // isRegionPoint
2950 pp,
2951 ppLocalPoints,
2952 snapDist,
2953 pointi,
2954 estimatedPt,
2955
2956 // Feature-point to pp point
2957 pointAttractor,
2959 // Feature-edge to pp point
2960 edgeAttractors,
2961 edgeConstraints,
2962 // pp point to nearest feature
2963 patchAttraction,
2964 patchConstraints
2965 );
2966 }
2967
2968 const pointIndexHit& info = nearInfo.second();
2969 if (featurePointStr && info.hit())
2970 {
2971 featurePointStr().writeLine(pt, info.point());
2972 }
2973 }
2974 }
2975 }
2976}
2977
2978
2979// Baffle handling
2980// ~~~~~~~~~~~~~~~
2981// Override pointAttractor, edgeAttractor, patchAttration etc. to
2982// implement 'baffle' handling.
2983// Baffle: the mesh pp point originates from a loose standing
2984// baffle.
2985// Sampling the surface with the surrounding face-centres only picks up
2986// a single triangle normal so above determineFeatures will not have
2987// detected anything. So explicitly pick up feature edges on the pp
2988// (after duplicating points & smoothing so will already have been
2989// expanded) and match these to the features.
2990void Foam::snappySnapDriver::determineBaffleFeatures
2991(
2992 const label iter,
2993 const bool baffleFeaturePoints,
2994 const scalar featureCos,
2995
2997 const pointField& ppLocalPoints,
2998 const scalarField& snapDist,
2999
3000 // Feature-point to pp point
3001 List<labelList>& pointAttractor,
3003 // Feature-edge to pp point
3004 List<List<DynamicList<point>>>& edgeAttractors,
3005 List<List<DynamicList<pointConstraint>>>& edgeConstraints,
3006 // pp point to nearest feature
3007 vectorField& patchAttraction,
3008 List<pointConstraint>& patchConstraints
3009) const
3010{
3011 const fvMesh& mesh = meshRefiner_.mesh();
3012 const refinementFeatures& features = meshRefiner_.features();
3013
3014 // Calculate edge-faces
3015 List<List<point>> edgeFaceNormals(pp.nEdges());
3016
3017 vectorField ppFaceNormals(pp.size());
3018 forAll(ppFaceNormals, facei)
3019 {
3020 ppFaceNormals[facei] = pp.localFaces()[facei].unitNormal(ppLocalPoints);
3021 }
3022
3023
3024 // Fill local data
3025 forAll(pp.edgeFaces(), edgei)
3026 {
3027 const labelList& eFaces = pp.edgeFaces()[edgei];
3028 List<point>& eFc = edgeFaceNormals[edgei];
3029 eFc.setSize(eFaces.size());
3030 forAll(eFaces, i)
3031 {
3032 label facei = eFaces[i];
3033 //eFc[i] = pp.faceNormals()[facei];
3034 eFc[i] = ppFaceNormals[facei];
3035 }
3036 }
3037
3038 {
3039 // Precalculate mesh edges for pp.edges.
3040 const labelList meshEdges
3041 (
3043 );
3044 // Collect all coupled edges. Does not filter duplicates/order
3046 (
3047 mesh,
3048 meshEdges,
3049 edgeFaceNormals,
3051 List<point>(),
3053 identityOp() // No flipping
3054 );
3055 }
3056
3057 // Detect baffle edges. Assume initial mesh will have 0,90 or 180
3058 // (baffle) degree angles so smoothing should make 0,90
3059 // to be less than 90. Choose reasonable value
3060 const scalar baffleFeatureCos = Foam::cos(degToRad(110.0));
3061
3062
3063 autoPtr<OBJstream> baffleEdgeStr;
3065 {
3066 baffleEdgeStr.reset
3067 (
3068 new OBJstream
3069 (
3070 meshRefiner_.mesh().time().path()
3071 / "baffleEdge_" + name(iter) + ".obj"
3072 )
3073 );
3074 Info<< nl << "Dumping baffle-edges to "
3075 << baffleEdgeStr().name() << endl;
3076 }
3077
3078
3079 // Is edge on baffle
3080 bitSet isBaffleEdge(pp.nEdges());
3081 label nBaffleEdges = 0;
3082 // Is point on
3083 // 0 : baffle-edge (0)
3084 // 1 : baffle-feature-point (1)
3085 // -1 : rest
3086 labelList pointStatus(pp.nPoints(), -1);
3087
3088 forAll(edgeFaceNormals, edgei)
3089 {
3090 const List<point>& efn = edgeFaceNormals[edgei];
3091
3092 if (efn.size() == 2 && (efn[0]&efn[1]) < baffleFeatureCos)
3093 {
3094 isBaffleEdge.set(edgei);
3095 ++nBaffleEdges;
3096 const edge& e = pp.edges()[edgei];
3097 pointStatus[e[0]] = 0;
3098 pointStatus[e[1]] = 0;
3099
3100 if (baffleEdgeStr)
3101 {
3102 baffleEdgeStr().write(e, ppLocalPoints);
3103 }
3104 }
3105 }
3106
3107 reduce(nBaffleEdges, sumOp<label>());
3108
3109 Info<< "Detected " << nBaffleEdges
3110 << " baffle edges out of "
3112 << " edges." << endl;
3113
3114
3115 //- Baffle edges will be too ragged to sensibly determine feature points
3116 //forAll(pp.pointEdges(), pointi)
3117 //{
3118 // if
3119 // (
3120 // isFeaturePoint
3121 // (
3122 // featureCos,
3123 // pp,
3124 // isBaffleEdge,
3125 // pointi
3126 // )
3127 // )
3128 // {
3129 // //Pout<< "Detected feature point:" << ppLocalPoints[pointi]
3130 // // << endl;
3131 // //-TEMPORARILY DISABLED:
3132 // //pointStatus[pointi] = 1;
3133 // }
3134 //}
3135
3136
3137 label nBafflePoints = 0;
3138 forAll(pointStatus, pointi)
3139 {
3140 if (pointStatus[pointi] != -1)
3141 {
3142 nBafflePoints++;
3143 }
3144 }
3145 reduce(nBafflePoints, sumOp<label>());
3146
3147
3148 label nPointAttract = 0;
3149 label nEdgeAttract = 0;
3150
3151 forAll(pointStatus, pointi)
3152 {
3153 const point& pt = ppLocalPoints[pointi];
3154
3155 if (pointStatus[pointi] == 0) // baffle edge
3156 {
3157 // 1: attract to near feature edge first
3158
3159 Tuple2<label, pointIndexHit> nearInfo = findNearFeatureEdge
3160 (
3161 false, // isRegionPoint?
3162 pp,
3163 ppLocalPoints,
3164 snapDist,
3165 pointi,
3166 pt,
3167
3168 edgeAttractors,
3169 edgeConstraints,
3170 patchAttraction,
3171 patchConstraints
3172 );
3173
3174
3175 //- MEJ:
3176 // 2: optionally override with nearest feature point.
3177 // On baffles we don't have enough normals to construct a feature
3178 // point so assume all feature edges are close to feature points
3179 if (nearInfo.second().hit())
3180 {
3181 nEdgeAttract++;
3182
3183 if (baffleFeaturePoints)
3184 {
3185 nearInfo = findNearFeaturePoint
3186 (
3187 false, // isRegionPoint,
3188
3189 pp,
3190 ppLocalPoints,
3191 snapDist,
3192 pointi,
3193 pt, // estimatedPt,
3194
3195 // Feature-point to pp point
3196 pointAttractor,
3198 // Feature-edge to pp point
3199 edgeAttractors,
3200 edgeConstraints,
3201 // pp point to nearest feature
3202 patchAttraction,
3203 patchConstraints
3204 );
3205
3206 if (nearInfo.first() != -1)
3207 {
3208 nEdgeAttract--;
3209 nPointAttract++;
3210 }
3211 }
3212 }
3213 }
3214 else if (pointStatus[pointi] == 1) // baffle point
3215 {
3216 labelList nearFeat;
3217 List<pointIndexHit> nearInfo;
3218 features.findNearestPoint
3219 (
3220 pointField(1, pt),
3221 scalarField(1, sqr(snapDist[pointi])),
3222 nearFeat,
3223 nearInfo
3224 );
3225
3226 label feati = nearFeat[0];
3227
3228 if (feati != -1)
3229 {
3230 nPointAttract++;
3231
3232 label featPointi = nearInfo[0].index();
3233 const point& featPt = nearInfo[0].hitPoint();
3234 scalar distSqr = featPt.distSqr(pt);
3235
3236 // Check if already attracted
3237 label oldPointi = pointAttractor[feati][featPointi];
3238
3239 if
3240 (
3241 oldPointi == -1
3242 || (
3243 distSqr
3244 < featPt.distSqr(ppLocalPoints[oldPointi])
3245 )
3246 )
3247 {
3248 pointAttractor[feati][featPointi] = pointi;
3249 pointConstraints[feati][featPointi].first() = 3;
3250 pointConstraints[feati][featPointi].second() = Zero;
3251
3252 // Store for later use
3253 patchAttraction[pointi] = featPt-pt;
3254 patchConstraints[pointi] =
3255 pointConstraints[feati][featPointi];
3256
3257 if (oldPointi != -1)
3258 {
3259 // The current point is closer so wins. Reset
3260 // the old point to attract to nearest edge
3261 // instead.
3262 findNearFeatureEdge
3263 (
3264 false, // isRegionPoint
3265 pp,
3266 ppLocalPoints,
3267 snapDist,
3268 oldPointi,
3269 ppLocalPoints[oldPointi],
3270
3271 edgeAttractors,
3272 edgeConstraints,
3273 patchAttraction,
3274 patchConstraints
3275 );
3276 }
3277 }
3278 else
3279 {
3280 // Make it fall through to check below
3281 feati = -1;
3282 }
3283 }
3284
3285 // Not found a feature point or another point is already
3286 // closer to that feature
3287 if (feati == -1)
3288 {
3289 //Pout<< "*** Falling back to finding nearest feature"
3290 // << " edge"
3291 // << " for baffle-feature-point " << pt
3292 // << endl;
3293
3294 Tuple2<label, pointIndexHit> nearInfo = findNearFeatureEdge
3295 (
3296 false, // isRegionPoint
3297 pp,
3298 ppLocalPoints,
3299 snapDist,
3300 pointi,
3301 pt, // starting point
3302
3303 edgeAttractors,
3304 edgeConstraints,
3305 patchAttraction,
3306 patchConstraints
3307 );
3308
3309 if (nearInfo.first() != -1)
3310 {
3311 nEdgeAttract++;
3312 }
3313 }
3314 }
3315 }
3316
3317 reduce(nPointAttract, sumOp<label>());
3318 reduce(nEdgeAttract, sumOp<label>());
3319
3320 Info<< "Baffle points : " << nBafflePoints
3321 << " of which attracted to :" << nl
3322 << " feature point : " << nPointAttract << nl
3323 << " feature edge : " << nEdgeAttract << nl
3324 << " rest : " << nBafflePoints-nPointAttract-nEdgeAttract
3325 << nl
3326 << endl;
3327}
3328
3329
3330void Foam::snappySnapDriver::reverseAttractMeshPoints
3331(
3332 const label iter,
3333
3335 const pointField& ppLocalPoints,
3336 const scalarField& snapDist,
3337
3338 // Feature-point to pp point
3339 const List<labelList>& pointAttractor,
3341 // Feature-edge to pp point
3342 const List<List<DynamicList<point>>>& edgeAttractors,
3343 const List<List<DynamicList<pointConstraint>>>& edgeConstraints,
3344
3345 const vectorField& rawPatchAttraction,
3346 const List<pointConstraint>& rawPatchConstraints,
3347
3348 // pp point to nearest feature
3349 vectorField& patchAttraction,
3350 List<pointConstraint>& patchConstraints
3351) const
3352{
3353 const refinementFeatures& features = meshRefiner_.features();
3354
3355 // Find nearest mesh point to feature edge
3356 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
3357 // Reverse lookup : go through all edgeAttractors and find the
3358 // nearest point on pp
3359
3360 // Get search domain and extend it a bit
3361 treeBoundBox bb(ppLocalPoints);
3362 {
3363 // Random number generator. Bit dodgy since not exactly random ;-)
3364 Random rndGen(65431);
3365
3366 // Slightly extended bb. Slightly off-centred just so on symmetric
3367 // geometry there are less face/edge aligned items.
3368 bb.inflate(rndGen, 1e-4, ROOTVSMALL);
3369 }
3370
3371 // Collect candidate points for attraction
3372 DynamicList<label> attractPoints(pp.nPoints());
3373 {
3374 const fvMesh& mesh = meshRefiner_.mesh();
3375
3376 boolList isFeatureEdgeOrPoint(pp.nPoints(), false);
3377 label nFeats = 0;
3378 forAll(rawPatchConstraints, pointi)
3379 {
3380 if (rawPatchConstraints[pointi].first() >= 2)
3381 {
3382 isFeatureEdgeOrPoint[pointi] = true;
3383 nFeats++;
3384 }
3385 }
3386
3387 Info<< "Initially selected " << returnReduce(nFeats, sumOp<label>())
3388 << " mesh points out of "
3390 << " for reverse attraction." << endl;
3391
3392 // Make sure is synchronised (note: check if constraint is already
3393 // synced in which case this is not needed here)
3395 (
3396 mesh,
3397 pp.meshPoints(),
3398 isFeatureEdgeOrPoint,
3399 orEqOp<bool>(), // combine op
3400 false
3401 );
3402
3403 for (label nGrow = 0; nGrow < 1; nGrow++)
3404 {
3405 boolList newIsFeatureEdgeOrPoint(isFeatureEdgeOrPoint);
3406
3407 forAll(pp.localFaces(), facei)
3408 {
3409 const face& f = pp.localFaces()[facei];
3410
3411 forAll(f, fp)
3412 {
3413 if (isFeatureEdgeOrPoint[f[fp]])
3414 {
3415 // Mark all points on face
3416 forAll(f, fp)
3417 {
3418 newIsFeatureEdgeOrPoint[f[fp]] = true;
3419 }
3420 break;
3421 }
3422 }
3423 }
3424
3425 isFeatureEdgeOrPoint = newIsFeatureEdgeOrPoint;
3426
3428 (
3429 mesh,
3430 pp.meshPoints(),
3431 isFeatureEdgeOrPoint,
3432 orEqOp<bool>(), // combine op
3433 false
3434 );
3435 }
3436
3437
3438 // Collect attractPoints
3439 forAll(isFeatureEdgeOrPoint, pointi)
3440 {
3441 if (isFeatureEdgeOrPoint[pointi])
3442 {
3443 attractPoints.append(pointi);
3444 }
3445 }
3446
3447 Info<< "Selected "
3448 << returnReduce(attractPoints.size(), sumOp<label>())
3449 << " mesh points out of "
3451 << " for reverse attraction." << endl;
3452 }
3453
3454
3456 (
3457 treeDataPoint(ppLocalPoints, attractPoints),
3458 bb, // overall search domain
3459 8, // maxLevel
3460 10, // leafsize
3461 3.0 // duplicity
3462 );
3463
3464 // Per mesh point the point on nearest feature edge.
3465 patchAttraction.setSize(pp.nPoints());
3466 patchAttraction = Zero;
3467 patchConstraints.setSize(pp.nPoints());
3468 patchConstraints = pointConstraint();
3469
3470 forAll(edgeAttractors, feati)
3471 {
3472 const List<DynamicList<point>>& edgeAttr = edgeAttractors[feati];
3473 const List<DynamicList<pointConstraint>>& edgeConstr =
3474 edgeConstraints[feati];
3475
3476 forAll(edgeAttr, featEdgei)
3477 {
3478 const DynamicList<point>& attr = edgeAttr[featEdgei];
3479 forAll(attr, i)
3480 {
3481 // Find nearest pp point
3482 const point& featPt = attr[i];
3483 pointIndexHit nearInfo = ppTree.findNearest
3484 (
3485 featPt,
3486 sqr(GREAT)
3487 );
3488
3489 if (nearInfo.hit())
3490 {
3491 label pointi =
3492 ppTree.shapes().objectIndex(nearInfo.index());
3493
3494 const point attraction
3495 (
3496 featPt
3497 - ppTree.shapes()[nearInfo.index()]
3498 );
3499
3500 // Check if this point is already being attracted. If so
3501 // override it only if nearer.
3502 if
3503 (
3504 patchConstraints[pointi].first() <= 1
3505 || magSqr(attraction) < magSqr(patchAttraction[pointi])
3506 )
3507 {
3508 patchAttraction[pointi] = attraction;
3509 patchConstraints[pointi] = edgeConstr[featEdgei][i];
3510 }
3511 }
3512 else
3513 {
3514 static label nWarn = 0;
3515
3516 if (nWarn < 100)
3517 {
3519 << "Did not find pp point near " << featPt
3520 << endl;
3521 nWarn++;
3522 if (nWarn == 100)
3523 {
3525 << "Reached warning limit " << nWarn
3526 << ". Suppressing further warnings." << endl;
3527 }
3528 }
3529
3530 }
3531 }
3532 }
3533 }
3534
3535
3536 // Different procs might have different patchAttraction,patchConstraints
3537 // however these only contain geometric information, no topology
3538 // so as long as we synchronise after overriding with feature points
3539 // there is no problem, just possibly a small error.
3540
3541
3542 // Find nearest mesh point to feature point
3543 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
3544 // (overrides attraction to feature edge)
3545 forAll(pointAttractor, feati)
3546 {
3547 const labelList& pointAttr = pointAttractor[feati];
3548 const List<pointConstraint>& pointConstr = pointConstraints[feati];
3549
3550 forAll(pointAttr, featPointi)
3551 {
3552 if (pointAttr[featPointi] != -1)
3553 {
3554 const point& featPt = features[feati].points()[featPointi];
3555
3556 // Find nearest pp point
3557 pointIndexHit nearInfo = ppTree.findNearest
3558 (
3559 featPt,
3560 sqr(GREAT)
3561 );
3562
3563 if (nearInfo.hit())
3564 {
3565 label pointi =
3566 ppTree.shapes().objectIndex(nearInfo.index());
3567
3568 const point attraction
3569 (
3570 featPt
3571 - ppTree.shapes()[nearInfo.index()]
3572 );
3573
3574 // - already attracted to feature edge : point always wins
3575 // - already attracted to feature point: nearest wins
3576
3577 if (patchConstraints[pointi].first() <= 1)
3578 {
3579 patchAttraction[pointi] = attraction;
3580 patchConstraints[pointi] = pointConstr[featPointi];
3581 }
3582 else if (patchConstraints[pointi].first() == 2)
3583 {
3584 patchAttraction[pointi] = attraction;
3585 patchConstraints[pointi] = pointConstr[featPointi];
3586 }
3587 else if (patchConstraints[pointi].first() == 3)
3588 {
3589 // Only if nearer
3590 if
3591 (
3592 magSqr(attraction)
3593 < magSqr(patchAttraction[pointi])
3594 )
3595 {
3596 patchAttraction[pointi] = attraction;
3597 patchConstraints[pointi] =
3598 pointConstr[featPointi];
3599 }
3600 }
3601 }
3602 }
3603 }
3604 }
3605}
3606
3607
3608void Foam::snappySnapDriver::featureAttractionUsingFeatureEdges
3609(
3610 const label iter,
3611 const bool multiRegionFeatureSnap,
3612 const bool strictRegionFeatureSnap, // special >=3 patch points handling
3613
3614 const bool detectBaffles,
3615 const bool baffleFeaturePoints,
3616
3617 const bool releasePoints,
3618 const bool stringFeatures,
3619 const bool avoidDiagonal,
3620
3621 const scalar featureCos,
3622
3624 const pointField& ppLocalPoints,
3625 const scalarField& snapDist,
3626 const vectorField& nearestDisp,
3627 const vectorField& nearestNormal,
3628
3629 const List<List<point>>& pointFaceSurfNormals,
3630 const List<List<point>>& pointFaceDisp,
3631 const List<List<point>>& pointFaceCentres,
3632 const labelListList& pointFacePatchID,
3633
3634 vectorField& patchAttraction,
3635 List<pointConstraint>& patchConstraints
3636) const
3637{
3638 const refinementFeatures& features = meshRefiner_.features();
3639 const fvMesh& mesh = meshRefiner_.mesh();
3640
3641 const bitSet isPatchMasterPoint
3642 (
3644 (
3645 mesh,
3646 pp.meshPoints()
3647 )
3648 );
3649
3650
3651 // Collect ordered attractions on feature edges
3652 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
3653
3654 // Per feature, per feature-edge a list of attraction points and their
3655 // originating vertex.
3656 List<List<DynamicList<point>>> edgeAttractors(features.size());
3658 (
3659 features.size()
3660 );
3661 forAll(features, feati)
3662 {
3663 label nFeatEdges = features[feati].edges().size();
3664 edgeAttractors[feati].setSize(nFeatEdges);
3665 edgeConstraints[feati].setSize(nFeatEdges);
3666 }
3667
3668 // Per feature, per feature-point the pp point that is attracted to it.
3669 // This list is only used to subset the feature-points that are actually
3670 // used.
3671 List<labelList> pointAttractor(features.size());
3673 forAll(features, feati)
3674 {
3675 label nFeatPoints = features[feati].points().size();
3676 pointAttractor[feati].setSize(nFeatPoints, -1);
3677 pointConstraints[feati].setSize(nFeatPoints);
3678 }
3679
3680 // Reverse: from pp point to nearest feature
3681 vectorField rawPatchAttraction(pp.nPoints(), Zero);
3682 List<pointConstraint> rawPatchConstraints(pp.nPoints());
3683
3684 determineFeatures
3685 (
3686 iter,
3687 featureCos,
3688 multiRegionFeatureSnap,
3689 strictRegionFeatureSnap, // extra logic for points on >=3 patches
3690
3691 pp,
3692 ppLocalPoints,
3693 snapDist, // per point max distance and nearest surface
3694 nearestDisp,
3695
3696 pointFaceSurfNormals, // per face nearest surface
3697 pointFaceDisp,
3698 pointFaceCentres,
3699 pointFacePatchID,
3700
3701 // Feature-point to pp point
3702 pointAttractor,
3704 // Feature-edge to pp point
3705 edgeAttractors,
3706 edgeConstraints,
3707 // pp point to nearest feature
3708 rawPatchAttraction,
3709 rawPatchConstraints
3710 );
3711
3712 // Print a bit about the attraction from patch point to feature
3713 if (debug)
3714 {
3715 Info<< "Raw geometric feature analysis : ";
3716 writeStats(pp, isPatchMasterPoint, rawPatchConstraints);
3717 }
3718
3719 // Baffle handling
3720 // ~~~~~~~~~~~~~~~
3721 // Override pointAttractor, edgeAttractor, rawPatchAttration etc. to
3722 // implement 'baffle' handling.
3723 // Baffle: the mesh pp point originates from a loose standing
3724 // baffle.
3725 // Sampling the surface with the surrounding face-centres only picks up
3726 // a single triangle normal so above determineFeatures will not have
3727 // detected anything. So explicitly pick up feature edges on the pp
3728 // (after duplicating points & smoothing so will already have been
3729 // expanded) and match these to the features.
3730 if (detectBaffles)
3731 {
3732 determineBaffleFeatures
3733 (
3734 iter,
3735 baffleFeaturePoints,
3736 featureCos,
3737
3738 pp,
3739 ppLocalPoints,
3740 snapDist,
3741
3742 // Feature-point to pp point
3743 pointAttractor,
3745 // Feature-edge to pp point
3746 edgeAttractors,
3747 edgeConstraints,
3748 // pp point to nearest feature
3749 rawPatchAttraction,
3750 rawPatchConstraints
3751 );
3752 }
3753
3754 // Print a bit about the attraction from patch point to feature
3755 if (debug)
3756 {
3757 Info<< "After baffle feature analysis : ";
3758 writeStats(pp, isPatchMasterPoint, rawPatchConstraints);
3759 }
3760
3761
3762 // Reverse lookup: Find nearest mesh point to feature edge
3763 // ~~~~~~~~~~~~~~~~----------------~~~~~~~~~~~~~~~~~~~~~~~
3764 // go through all edgeAttractors and find the nearest point on pp
3765
3766 reverseAttractMeshPoints
3767 (
3768 iter,
3769
3770 pp,
3771 ppLocalPoints,
3772 snapDist,
3773
3774 // Feature-point to pp point
3775 pointAttractor,
3777 // Feature-edge to pp point
3778 edgeAttractors,
3779 edgeConstraints,
3780
3781 // Estimated feature point
3782 rawPatchAttraction,
3783 rawPatchConstraints,
3784
3785 // pp point to nearest feature
3786 patchAttraction,
3787 patchConstraints
3788 );
3789
3790 // Print a bit about the attraction from patch point to feature
3791 if (debug)
3792 {
3793 Info<< "Reverse attract feature analysis : ";
3794 writeStats(pp, isPatchMasterPoint, patchConstraints);
3795 }
3796
3797 // Dump
3799 {
3800 OBJstream featureEdgeStr
3801 (
3802 meshRefiner_.mesh().time().path()
3803 / "edgeAttractors_" + name(iter) + ".obj"
3804 );
3805 Info<< "Dumping feature-edge attraction to "
3806 << featureEdgeStr.name() << endl;
3807
3808 OBJstream featurePointStr
3809 (
3810 meshRefiner_.mesh().time().path()
3811 / "pointAttractors_" + name(iter) + ".obj"
3812 );
3813 Info<< "Dumping feature-point attraction to "
3814 << featurePointStr.name() << endl;
3815
3816 forAll(patchConstraints, pointi)
3817 {
3818 const point& pt = ppLocalPoints[pointi];
3819 const vector& attr = patchAttraction[pointi];
3820
3821 if (patchConstraints[pointi].first() == 2)
3822 {
3823 featureEdgeStr.writeLine(pt, pt+attr);
3824 }
3825 else if (patchConstraints[pointi].first() == 3)
3826 {
3827 featurePointStr.writeLine(pt, pt+attr);
3828 }
3829 }
3830 }
3831
3832
3833 //MEJ: any faces that have multi-patch points only keep the
3834 // multi-patch
3835 // points. The other points on the face will be dragged along
3836 // (hopefully)
3837 if (releasePoints)
3838 {
3839 releasePointsNextToMultiPatch
3840 (
3841 iter,
3842 featureCos,
3843
3844 pp,
3845 ppLocalPoints,
3846 snapDist,
3847
3848 pointFaceCentres,
3849 pointFacePatchID,
3850
3851 rawPatchAttraction,
3852 rawPatchConstraints,
3853
3854 patchAttraction,
3855 patchConstraints
3856 );
3857 }
3858
3859
3860 // Snap edges to feature edges
3861 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
3862 // Walk existing edges and snap remaining ones (that are marked as
3863 // feature edges in rawPatchConstraints)
3864 if (stringFeatures)
3865 {
3866 stringFeatureEdges
3867 (
3868 iter,
3869 featureCos,
3870
3871 pp,
3872 ppLocalPoints,
3873 snapDist,
3874
3875 rawPatchAttraction,
3876 rawPatchConstraints,
3877
3878 patchAttraction,
3879 patchConstraints
3880 );
3881 }
3882
3883
3884 // Avoid diagonal attraction
3885 // ~~~~~~~~~~~~~~~~~~~~~~~~~
3886 // Attract one of the non-diagonal points.
3887 if (avoidDiagonal)
3888 {
3889 avoidDiagonalAttraction
3890 (
3891 iter,
3892 featureCos,
3893 pp,
3894 ppLocalPoints,
3895 patchAttraction,
3896 patchConstraints
3897 );
3898 }
3899
3900
3902 {
3903 dumpMove
3904 (
3905 meshRefiner_.mesh().time().path()
3906 / "patchAttraction_" + name(iter) + ".obj",
3907 ppLocalPoints,
3908 ppLocalPoints + patchAttraction
3909 );
3910 }
3911}
3912
3913
3914// Correct for squeezing of face
3915void Foam::snappySnapDriver::preventFaceSqueeze
3916(
3917 const label iter,
3918 const scalar featureCos,
3919
3921 const pointField& ppLocalPoints,
3922 const scalarField& snapDist,
3923 const vectorField& nearestAttraction,
3924
3925 vectorField& patchAttraction,
3926 List<pointConstraint>& patchConstraints
3927) const
3928{
3929 autoPtr<OBJstream> strPtr;
3931 {
3932 strPtr.reset
3933 (
3934 new OBJstream
3935 (
3936 meshRefiner_.mesh().time().path()
3937 / "faceSqueeze_" + name(iter) + ".obj"
3938 )
3939 );
3940 Info<< "Dumping faceSqueeze corrections to "
3941 << strPtr().name() << endl;
3942 }
3943
3945 face singleF;
3946 forAll(pp.localFaces(), facei)
3947 {
3948 const face& f = pp.localFaces()[facei];
3949
3950 if (f.size() != points.size())
3951 {
3952 points.setSize(f.size());
3953 singleF.setSize(f.size());
3954 for (label i = 0; i < f.size(); i++)
3955 {
3956 singleF[i] = i;
3957 }
3958 }
3959 label nConstraints = 0;
3960 forAll(f, fp)
3961 {
3962 label pointi = f[fp];
3963 const point& pt = ppLocalPoints[pointi];
3964
3965 if (patchConstraints[pointi].first() > 1)
3966 {
3967 points[fp] = pt + patchAttraction[pointi];
3968 nConstraints++;
3969 }
3970 else
3971 {
3972 points[fp] = pt;
3973 }
3974 }
3975
3976 if (nConstraints == f.size())
3977 {
3978 if (f.size() == 3)
3979 {
3980 // Triangle: knock out attraction altogether
3981
3982 // For now keep the points on the longest edge
3983 label maxFp = -1;
3984 scalar maxS = -1;
3985 forAll(f, fp)
3986 {
3987 const point& pt = ppLocalPoints[f[fp]];
3988 const point& nextPt = ppLocalPoints[f.nextLabel(fp)];
3989
3990 scalar s = pt.distSqr(nextPt);
3991 if (s > maxS)
3992 {
3993 maxS = s;
3994 maxFp = fp;
3995 }
3996 }
3997 if (maxFp != -1)
3998 {
3999 label pointi = f.prevLabel(maxFp);
4000
4001 // Reset attraction on pointi to nearest
4002
4003 const point& pt = ppLocalPoints[pointi];
4004
4005 //Pout<< "** on triangle " << pp.faceCentres()[facei]
4006 // << " knocking out attraction to " << pointi
4007 // << " at:" << pt
4008 // << endl;
4009
4010 patchAttraction[pointi] = nearestAttraction[pointi];
4011
4012 if (strPtr)
4013 {
4014 strPtr().writeLine(pt, pt+patchAttraction[pointi]);
4015 }
4016 }
4017 }
4018 else
4019 {
4020 scalar oldArea = f.mag(ppLocalPoints);
4021 scalar newArea = singleF.mag(points);
4022 if (newArea < 0.1*oldArea)
4023 {
4024 // For now remove the point with largest distance
4025 label maxFp = -1;
4026 scalar maxS = -1;
4027 forAll(f, fp)
4028 {
4029 scalar s = magSqr(patchAttraction[f[fp]]);
4030 if (s > maxS)
4031 {
4032 maxS = s;
4033 maxFp = fp;
4034 }
4035 }
4036 if (maxFp != -1)
4037 {
4038 label pointi = f[maxFp];
4039 // Lower attraction on pointi
4040 patchAttraction[pointi] *= 0.5;
4041 }
4042 }
4043 }
4044 }
4045 }
4046}
4047
4048
4049Foam::vectorField Foam::snappySnapDriver::calcNearestSurfaceFeature
4050(
4051 const snapParameters& snapParams,
4052 const bool alignMeshEdges,
4053 const bool strictRegionFeatureSnap, // use patches
4054 const label iter,
4055 const scalar featureCos,
4056 const scalar featureAttract,
4057 const scalarField& snapDist,
4058 const vectorField& nearestDisp,
4059 const vectorField& nearestNormal,
4061 const pointField& ppLocalPoints,
4062 vectorField& patchAttraction,
4063 List<pointConstraint>& patchConstraints,
4064
4065 DynamicList<label>& splitFaces,
4066 DynamicList<labelPair>& splits,
4067 DynamicList<labelPair>& splitPatches
4068) const
4069{
4070 if (dryRun_)
4071 {
4072 return nearestDisp;
4073 }
4074
4075
4076 const Switch implicitFeatureAttraction = snapParams.implicitFeatureSnap();
4077 const Switch explicitFeatureAttraction = snapParams.explicitFeatureSnap();
4078 const Switch multiRegionFeatureSnap = snapParams.multiRegionFeatureSnap();
4079
4080 Info<< "Overriding displacement on features :" << nl
4081 << " implicit features : " << implicitFeatureAttraction << nl
4082 << " explicit features : " << explicitFeatureAttraction << nl
4083 << " multi-patch features : " << multiRegionFeatureSnap << nl
4084 << endl;
4085
4086 const fvMesh& mesh = meshRefiner_.mesh();
4087
4088
4089 //const bitSet isMasterPoint(syncTools::getMasterPoints(mesh));
4090 const bitSet isPatchMasterPoint
4091 (
4093 (
4094 mesh,
4095 pp.meshPoints()
4096 )
4097 );
4098
4099 // Per point, per surrounding face:
4100 // - faceSurfaceNormal
4101 // - faceDisp
4102 // - faceCentres
4103 List<List<point>> pointFaceSurfNormals;
4104 List<List<point>> pointFaceDisp;
4105 List<List<point>> pointFaceCentres;
4106 List<labelList> pointFacePatchID;
4107
4108 {
4109 // Calculate attraction distance per face (from the attraction distance
4110 // per point)
4111 scalarField faceSnapDist(pp.size(), -GREAT);
4112 forAll(pp.localFaces(), facei)
4113 {
4114 const face& f = pp.localFaces()[facei];
4115 forAll(f, fp)
4116 {
4117 faceSnapDist[facei] = max(faceSnapDist[facei], snapDist[f[fp]]);
4118 }
4119 }
4120
4121
4122 // Displacement and orientation per pp face
4123 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
4124
4125 // vector from point on surface back to face centre
4126 vectorField faceDisp(pp.size(), Zero);
4127 // normal of surface at point on surface
4128 vectorField faceSurfaceNormal(pp.size(), Zero);
4129 labelList faceSurfaceGlobalRegion(pp.size(), -1);
4130 //vectorField faceRotation(pp.size(), Zero);
4131
4132 calcNearestFace
4133 (
4134 iter,
4135 pp,
4136 ppLocalPoints,
4137 faceSnapDist,
4138 faceDisp,
4139 faceSurfaceNormal,
4140 faceSurfaceGlobalRegion
4141 //faceRotation
4142 );
4143
4144
4145 // Collect (possibly remote) per point data of all surrounding faces
4146 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
4147 // - faceSurfaceNormal
4148 // - faceDisp
4149 // - faceCentres
4150 calcNearestFacePointProperties
4151 (
4152 iter,
4153 pp,
4154 ppLocalPoints,
4155 faceDisp,
4156 faceSurfaceNormal,
4157 faceSurfaceGlobalRegion,
4158
4159 pointFaceSurfNormals,
4160 pointFaceDisp,
4161 pointFaceCentres,
4162 pointFacePatchID
4163 );
4164 }
4165
4166
4167 // Start off with nearest point on surface
4168 vectorField patchDisp = nearestDisp;
4169
4170
4171 // Main calculation
4172 // ~~~~~~~~~~~~~~~~
4173 // This is the main intelligence which calculates per point the vector to
4174 // attract it to the nearest surface. There are lots of possibilities
4175 // here.
4176
4177 // Nearest feature
4178 patchAttraction.setSize(ppLocalPoints.size());
4179 patchAttraction = Zero;
4180 // Constraints at feature
4181 patchConstraints.setSize(ppLocalPoints.size());
4182 patchConstraints = pointConstraint();
4183
4184 if (implicitFeatureAttraction)
4185 {
4186 // Sample faces around each point and see if nearest surface normal
4187 // differs. Reconstruct a feature edge/point if possible and snap to
4188 // it.
4189 featureAttractionUsingReconstruction
4190 (
4191 iter,
4192 featureCos,
4193
4194 pp,
4195 ppLocalPoints,
4196 snapDist,
4197 nearestDisp,
4198
4199 pointFaceSurfNormals,
4200 pointFaceDisp,
4201 pointFaceCentres,
4202 pointFacePatchID,
4203
4204 patchAttraction,
4205 patchConstraints
4206 );
4207 }
4208
4209 if (explicitFeatureAttraction)
4210 {
4211 // Only do fancy stuff if alignMeshEdges
4212 bool releasePoints = false;
4213 bool stringFeatures = false;
4214 bool avoidDiagonal = false;
4215 if (alignMeshEdges)
4216 {
4217 releasePoints = snapParams.releasePoints();
4218 stringFeatures = snapParams.stringFeatures();
4219 avoidDiagonal = snapParams.avoidDiagonal();
4220 }
4221
4222
4223 // Sample faces around each point and see if nearest surface normal
4224 // differs. For those find the nearest real feature edge/point and
4225 // store the correspondence. Then loop over feature edge/point
4226 // and attract those nearest mesh point. (the first phase just is
4227 // a subsetting of candidate points, the second makes sure that only
4228 // one mesh point gets attracted per feature)
4229 featureAttractionUsingFeatureEdges
4230 (
4231 iter,
4232 multiRegionFeatureSnap,
4233 strictRegionFeatureSnap, // special logic for >=3 patches points?
4234
4235 snapParams.detectBaffles(),
4236 snapParams.baffleFeaturePoints(), // all points on baffle edges
4237 // are attracted to feature pts
4238
4239 releasePoints,
4240 stringFeatures,
4241 avoidDiagonal,
4242
4243 featureCos,
4244
4245 pp,
4246 ppLocalPoints,
4247 snapDist,
4248 nearestDisp,
4249 nearestNormal,
4250
4251 pointFaceSurfNormals,
4252 pointFaceDisp,
4253 pointFaceCentres,
4254 pointFacePatchID,
4255
4256 patchAttraction,
4257 patchConstraints
4258 );
4259 }
4260
4261 if (!alignMeshEdges)
4262 {
4263 const scalar concaveCos = Foam::cos
4264 (
4265 degToRad(snapParams.concaveAngle())
4266 );
4267 const scalar minAreaRatio = snapParams.minAreaRatio();
4268
4269 Info<< "Introducing face splits to avoid rotating"
4270 << " mesh edges. Splitting faces when" << nl
4271 << indent << "- angle not concave by more than "
4272 << snapParams.concaveAngle() << " degrees" << nl
4273 << indent << "- resulting triangles of similar area "
4274 << " (ratio within " << minAreaRatio << ")" << nl
4275 << endl;
4276
4277 splitDiagonals
4278 (
4279 featureCos,
4280 concaveCos,
4281 minAreaRatio,
4282 pp,
4283 ppLocalPoints,
4284
4285 nearestDisp,
4286 nearestNormal,
4287 pointFacePatchID, // per point the connected patches. Is used
4288 // to re-patch newly split faces
4289 patchAttraction,
4290 patchConstraints,
4291
4292 splitFaces,
4293 splits,
4294 splitPatches
4295 );
4296
4297 if (debug)
4298 {
4299 Info<< "Diagonal attraction feature correction : ";
4300 writeStats(pp, isPatchMasterPoint, patchConstraints);
4301 }
4302 }
4303
4304
4305 preventFaceSqueeze
4306 (
4307 iter,
4308 featureCos,
4309
4310 pp,
4311 ppLocalPoints,
4312 snapDist,
4313 nearestDisp,
4314
4315 patchAttraction,
4316 patchConstraints
4317 );
4318
4319 {
4320 vector avgPatchDisp = meshRefinement::gAverage
4321 (
4322 isPatchMasterPoint,
4323 patchDisp
4324 );
4325 vector avgPatchAttr = meshRefinement::gAverage
4326 (
4327 isPatchMasterPoint,
4328 patchAttraction
4329 );
4330
4331 Info<< "Attraction:" << endl
4332 << " linear : max:" << gMaxMagSqr(patchDisp)
4333 << " avg:" << avgPatchDisp << endl
4334 << " feature : max:" << gMaxMagSqr(patchAttraction)
4335 << " avg:" << avgPatchAttr << endl;
4336 }
4337
4338 // So now we have:
4339 // - patchDisp : point movement to go to nearest point on surface
4340 // (either direct or through interpolation of
4341 // face nearest)
4342 // - patchAttraction : direct attraction to features
4343 // - patchConstraints : type of features
4344
4345 // Use any combination of patchDisp and direct feature attraction.
4346
4347
4348 // Mix with direct feature attraction
4349 forAll(patchConstraints, pointi)
4350 {
4351 if (patchConstraints[pointi].first() > 1)
4352 {
4353 patchDisp[pointi] =
4354 (1.0-featureAttract)*patchDisp[pointi]
4355 + featureAttract*patchAttraction[pointi];
4356 }
4357 }
4358
4359
4360
4361 // Count
4362 {
4363 Info<< "Feature analysis : ";
4364 writeStats(pp, isPatchMasterPoint, patchConstraints);
4365 }
4366
4367
4368 // Now we have the displacement per patch point to move onto the surface
4369 // Split into tangential and normal direction.
4370 // - start off with all non-constrained points following the constrained
4371 // ones since point normals not relevant.
4372 // - finish with only tangential component smoothed.
4373 // Note: tangential is most
4374 // likely to come purely from face-centre snapping, not face rotation.
4375 // Note: could use the constraints here (constraintTransformation())
4376 // but this is not necessarily accurate and we're smoothing to
4377 // get out of problems.
4378
4379 if (featureAttract < 1-0.001)
4380 {
4381 //const bitSet isMasterEdge(syncTools::getMasterEdges(mesh));
4382 const labelList meshEdges
4383 (
4385 );
4386 const bitSet isPatchMasterEdge
4387 (
4389 (
4390 mesh,
4391 meshEdges
4392 )
4393 );
4394
4395 const vectorField pointNormals
4396 (
4398 (
4399 mesh,
4400 pp
4401 )
4402 );
4403
4404 // 1. Smoothed all displacement
4405 vectorField smoothedPatchDisp = patchDisp;
4406 smoothAndConstrain
4407 (
4408 isPatchMasterEdge,
4409 pp,
4410 meshEdges,
4411 patchConstraints,
4412 smoothedPatchDisp
4413 );
4414
4415
4416 // 2. Smoothed tangential component
4417 vectorField tangPatchDisp = patchDisp;
4418 tangPatchDisp -= (pointNormals & patchDisp) * pointNormals;
4419 smoothAndConstrain
4420 (
4421 isPatchMasterEdge,
4422 pp,
4423 meshEdges,
4424 patchConstraints,
4425 tangPatchDisp
4426 );
4427
4428 // Re-add normal component
4429 tangPatchDisp += (pointNormals & patchDisp) * pointNormals;
4430
4432 {
4433 dumpMove
4434 (
4435 mesh.time().path()
4436 / "tangPatchDispConstrained_" + name(iter) + ".obj",
4437 ppLocalPoints,
4438 ppLocalPoints + tangPatchDisp
4439 );
4440 }
4441
4442 patchDisp =
4443 (1.0-featureAttract)*smoothedPatchDisp
4444 + featureAttract*tangPatchDisp;
4445 }
4446
4447
4448 const scalar relax = featureAttract;
4449 patchDisp *= relax;
4450
4451
4452 // Points on zones in one domain but only present as point on other
4453 // will not do condition 2 on all. Sync explicitly.
4455 (
4456 mesh,
4457 pp.meshPoints(),
4458 patchDisp,
4459 minMagSqrEqOp<point>(), // combine op
4460 vector(GREAT, GREAT, GREAT) // null value (note: cannot use VGREAT)
4461 );
4462
4463 return patchDisp;
4464}
4465
4466
4467// ************************************************************************* //
scalar y
constexpr scalar pi(M_PI)
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
const polyBoundaryMesh & pbm
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
const Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in the list.
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
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const bitSet &flipMap=bitSet::null())
Return parallel consistent point normals for patches using mesh points.
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const labelListList & pointFaces() const
Return point-face addressing.
const labelListList & edgeFaces() const
Return edge-face addressing.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Random number generator.
Definition Random.H:56
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
T & first()
Access first element of the list, position [0].
Definition UList.H:957
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form zero
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Non-pointer based hierarchical recursive searching.
void operator()(List< T > &x, const List< T > &y) const
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
Default transformation behaviour.
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
static void splitFace(const face &f, const labelPair &split, face &f0, face &f1)
Helper: split face into:
static bitSet getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
const fvMesh & mesh() const
Reference to mesh.
A reference point and direction.
Definition plane.H:114
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition plane.H:91
Accumulates point constraints through successive applications of the applyConstraint function.
Application of (multi-)patch point constraints.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const labelList & patchID() const
Per boundary face label the patch index.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
label nInternalFaces() const noexcept
Number of internal faces.
label nPoints() const noexcept
Number of mesh points.
label nFaces() const noexcept
Number of mesh faces.
const vectorField & faceAreas() const
Encapsulates queries for features.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
const labelList & surfaces() const
Simple container to keep together snap specific information.
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName).
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of named surfaces (surfaces with faceZoneName).
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition syncTools.C:119
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
Standard boundBox with extra functionality for use in octree.
Holds (reference to) pointField. Encapsulation of data needed for octree searches....
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
UEqn relax()
const volScalarField & p0
Definition EEqn.H:36
const polyBoundaryMesh & patches
static bool split(const std::string &line, std::string &key, std::string &val)
Definition cpuInfo.C:32
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
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.
const char * end
Definition SVGTools.H:223
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
void writeStats(Ostream &os, const extendedFeatureEdgeMesh &emesh)
Write some information.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
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.
Type gMaxMagSqr(const UList< Type > &f, const label comm)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
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
A functor that returns its argument unchanged (cf. C++20 std::identity) Should never be specialized.
Definition stdFoam.H:108
Unit conversion functions.
Random rndGen
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))