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