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