surfaceHookUp.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  surfaceHookUp
26 
27 Description
28  Find close open edges and stitches the surface along them
29 
30 Usage
31  - surfaceHookUp hookDistance [OPTION]
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "argList.H"
36 #include "Time.H"
37 
38 #include "triSurfaceMesh.H"
39 #include "indexedOctree.H"
40 #include "treeBoundBox.H"
41 #include "PackedBoolList.H"
42 #include "unitConversion.H"
43 #include "searchableSurfaces.H"
44 #include "systemDict.H"
45 
46 using namespace Foam;
47 
48 // Split facei along edgeI at position newPointi
49 void greenRefine
50 (
51  const triSurface& surf,
52  const label facei,
53  const label edgeI,
54  const label newPointi,
55  DynamicList<labelledTri>& newFaces
56 )
57 {
58  const labelledTri& f = surf.localFaces()[facei];
59  const edge& e = surf.edges()[edgeI];
60 
61  // Find index of edge in face.
62 
63  label fp0 = findIndex(f, e[0]);
64  label fp1 = f.fcIndex(fp0);
65  label fp2 = f.fcIndex(fp1);
66 
67  if (f[fp1] == e[1])
68  {
69  // Edge oriented like face
70  newFaces.append
71  (
73  (
74  f[fp0],
75  newPointi,
76  f[fp2],
77  f.region()
78  )
79  );
80  newFaces.append
81  (
83  (
84  newPointi,
85  f[fp1],
86  f[fp2],
87  f.region()
88  )
89  );
90  }
91  else
92  {
93  newFaces.append
94  (
96  (
97  f[fp2],
98  newPointi,
99  f[fp1],
100  f.region()
101  )
102  );
103  newFaces.append
104  (
106  (
107  newPointi,
108  f[fp0],
109  f[fp1],
110  f.region()
111  )
112  );
113  }
114 }
115 
116 
117 //scalar checkEdgeAngle
118 //(
119 // const triSurface& surf,
120 // const label edgeIndex,
121 // const label pointIndex,
122 // const scalar& angle
123 //)
124 //{
125 // const edge& e = surf.edges()[edgeIndex];
126 
127 // vector eVec = e.vec(surf.localPoints());
128 // eVec /= mag(eVec) + small;
129 
130 // const labelList& pEdges = surf.pointEdges()[pointIndex];
131 //
132 // forAll(pEdges, eI)
133 // {
134 // const edge& nearE = surf.edges()[pEdges[eI]];
135 
136 // vector nearEVec = nearE.vec(surf.localPoints());
137 // nearEVec /= mag(nearEVec) + small;
138 
139 // const scalar dot = eVec & nearEVec;
140 // const scalar minCos = degToRad(angle);
141 
142 // if (mag(dot) > minCos)
143 // {
144 // return false;
145 // }
146 // }
147 
148 // return true;
149 //}
150 
151 
152 void createBoundaryEdgeTrees
153 (
154  const PtrList<triSurfaceMesh>& surfs,
156  labelListList& treeBoundaryEdges
157 )
158 {
159  forAll(surfs, surfI)
160  {
161  const triSurface& surf = surfs[surfI];
162 
163  // Boundary edges
164  treeBoundaryEdges[surfI] =
165  labelList
166  (
167  identityMap(surf.nEdges() - surf.nInternalEdges())
168  + surf.nInternalEdges()
169  );
170 
171  Random rndGen(17301893);
172 
173  // Slightly extended bb. Slightly off-centred just so on symmetric
174  // geometry there are less face/edge aligned items.
175  treeBoundBox bb
176  (
178  );
179 
180  bEdgeTrees.set
181  (
182  surfI,
184  (
186  (
187  false, // cachebb
188  surf.edges(), // edges
189  surf.localPoints(), // points
190  treeBoundaryEdges[surfI] // selected edges
191  ),
192  bb, // bb
193  8, // maxLevel
194  10, // leafsize
195  3.0 // duplicity
196  )
197  );
198  }
199 }
200 
201 
202 class findNearestOpSubset
203 {
204  const indexedOctree<treeDataEdge>& tree_;
205 
206  DynamicList<label>& shapeMask_;
207 
208 public:
209 
210  findNearestOpSubset
211  (
212  const indexedOctree<treeDataEdge>& tree,
213  DynamicList<label>& shapeMask
214  )
215  :
216  tree_(tree),
217  shapeMask_(shapeMask)
218  {}
219 
220  void operator()
221  (
222  const labelUList& indices,
223  const point& sample,
224 
225  scalar& nearestDistSqr,
226  label& minIndex,
227  point& nearestPoint
228  ) const
229  {
230  const treeDataEdge& shape = tree_.shapes();
231 
232  forAll(indices, i)
233  {
234  const label index = indices[i];
235  const label edgeIndex = shape.edgeLabels()[index];
236 
237  if
238  (
239  !shapeMask_.empty()
240  && findIndex(shapeMask_, edgeIndex) != -1
241  )
242  {
243  continue;
244  }
245 
246  const edge& e = shape.edges()[edgeIndex];
247 
248  pointHit nearHit = e.line(shape.points()).nearestDist(sample);
249 
250  // Only register hit if closest point is not an edge point
251  if (nearHit.hit())
252  {
253  scalar distSqr = sqr(nearHit.distance());
254 
255  if (distSqr < nearestDistSqr)
256  {
257  nearestDistSqr = distSqr;
258  minIndex = index;
259  nearestPoint = nearHit.rawPoint();
260  }
261  }
262  }
263  }
264 };
265 
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 int main(int argc, char *argv[])
270 {
272  (
273  "hook surfaces to other surfaces by moving and retriangulating their"
274  "boundary edges to match other surface boundary edges"
275  );
277  argList::validArgs.append("hookTolerance");
278 
279  #include "addDictOption.H"
280 
281  #include "setRootCase.H"
282  #include "createTime.H"
283 
284  const dictionary dict(systemDict("surfaceHookUpDict", args, runTime));
285 
286  const scalar dist(args.argRead<scalar>(1));
287  const scalar matchTolerance(max(1e-6*dist, small));
288  const label maxIters = 100;
289 
290  Info<< "Hooking distance = " << dist << endl;
291 
292  searchableSurfaces surfs
293  (
294  IOobject
295  (
296  "surfacesToHook",
297  runTime.constant(),
299  runTime
300  ),
301  dict,
302  true // assume single-region names get surface name
303  );
304 
305  Info<< nl << "Reading surfaces: " << endl;
306  forAll(surfs, surfI)
307  {
308  Info<< incrIndent;
309  Info<< nl << indent << "Surface = " << surfs.names()[surfI] << endl;
310 
311  const wordList& regions = surfs[surfI].regions();
312  forAll(regions, surfRegionI)
313  {
314  Info<< incrIndent;
315  Info<< indent << "Regions = " << regions[surfRegionI] << endl;
316  Info<< decrIndent;
317  }
318  Info<< decrIndent;
319  }
320 
321  PtrList<indexedOctree<treeDataEdge>> bEdgeTrees(surfs.size());
322  labelListList treeBoundaryEdges(surfs.size());
323 
324  List<DynamicList<labelledTri>> newFaces(surfs.size());
325  List<DynamicList<point>> newPoints(surfs.size());
326  List<PackedBoolList> visitedFace(surfs.size());
327 
328  PtrList<triSurfaceMesh> newSurfaces(surfs.size());
329  forAll(surfs, surfI)
330  {
331  const triSurfaceMesh& surf =
332  refCast<const triSurfaceMesh>(surfs[surfI]);
333 
334  newSurfaces.set
335  (
336  surfI,
337  new triSurfaceMesh
338  (
339  IOobject
340  (
341  "hookedSurface_" + surfs.names()[surfI],
342  runTime.constant(),
344  runTime
345  ),
346  surf
347  )
348  );
349  }
350 
351  label nChanged = 0;
352  label nIters = 1;
353 
354  do
355  {
356  Info<< nl << "Iteration = " << nIters++ << endl;
357  nChanged = 0;
358 
359  createBoundaryEdgeTrees(newSurfaces, bEdgeTrees, treeBoundaryEdges);
360 
361  forAll(newSurfaces, surfI)
362  {
363  const triSurface& newSurf = newSurfaces[surfI];
364 
365  newFaces[surfI] = newSurf.localFaces();
366  newPoints[surfI] = newSurf.localPoints();
367  visitedFace[surfI] = PackedBoolList(newSurf.size(), false);
368  }
369 
370  forAll(newSurfaces, surfI)
371  {
372  const triSurface& surf = newSurfaces[surfI];
373 
374  List<pointIndexHit> bPointsTobEdges(surf.boundaryPoints().size());
375  labelList bPointsHitTree(surf.boundaryPoints().size(), -1);
376 
377  const labelListList& pointEdges = surf.pointEdges();
378 
379  forAll(bPointsTobEdges, bPointi)
380  {
381  pointIndexHit& nearestHit = bPointsTobEdges[bPointi];
382 
383  const label pointi = surf.boundaryPoints()[bPointi];
384  const point& samplePt = surf.localPoints()[pointi];
385 
386  const labelList& pEdges = pointEdges[pointi];
387 
388  // Add edges connected to the edge to the shapeMask
389  DynamicList<label> shapeMask;
390  shapeMask.append(pEdges);
391 
392  forAll(bEdgeTrees, treeI)
393  {
394  const indexedOctree<treeDataEdge>& bEdgeTree =
395  bEdgeTrees[treeI];
396 
397  pointIndexHit currentHit =
398  bEdgeTree.findNearest
399  (
400  samplePt,
401  sqr(dist),
402  findNearestOpSubset
403  (
404  bEdgeTree,
405  shapeMask
406  )
407  );
408 
409  if
410  (
411  currentHit.hit()
412  &&
413  (
414  !nearestHit.hit()
415  ||
416  (
417  magSqr(currentHit.hitPoint() - samplePt)
418  < magSqr(nearestHit.hitPoint() - samplePt)
419  )
420  )
421  )
422  {
423  nearestHit = currentHit;
424  bPointsHitTree[bPointi] = treeI;
425  }
426  }
427 
428  scalar dist2 = magSqr(nearestHit.rawPoint() - samplePt);
429 
430  if (nearestHit.hit())
431  {
432  // bool rejectEdge =
433  // checkEdgeAngle
434  // (
435  // surf,
436  // nearestHit.index(),
437  // pointi,
438  // 30
439  // );
440 
441  if (dist2 > Foam::sqr(dist))
442  {
443  nearestHit.setMiss();
444  }
445  }
446  }
447 
448  forAll(bPointsTobEdges, bPointi)
449  {
450  const pointIndexHit& eHit = bPointsTobEdges[bPointi];
451 
452  if (eHit.hit())
453  {
454  const label hitSurfI = bPointsHitTree[bPointi];
455  const triSurface& hitSurf = newSurfaces[hitSurfI];
456 
457  const label eIndex =
458  treeBoundaryEdges[hitSurfI][eHit.index()];
459  const edge& e = hitSurf.edges()[eIndex];
460 
461  const label pointi = surf.boundaryPoints()[bPointi];
462 
463  const labelList& eFaces = hitSurf.edgeFaces()[eIndex];
464 
465  if (eFaces.size() != 1)
466  {
468  << "Edge is attached to " << eFaces.size()
469  << " faces." << endl;
470 
471  continue;
472  }
473 
474  const label facei = eFaces[0];
475 
476  if (visitedFace[hitSurfI][facei])
477  {
478  continue;
479  }
480 
481  DynamicList<labelledTri> newFacesFromSplit(2);
482 
483  const point& pt = surf.localPoints()[pointi];
484 
485  if
486  (
487  (
488  magSqr(pt - hitSurf.localPoints()[e.start()])
489  < matchTolerance
490  )
491  || (
492  magSqr(pt - hitSurf.localPoints()[e.end()])
493  < matchTolerance
494  )
495  )
496  {
497  continue;
498  }
499 
500  nChanged++;
501 
502  label newPointi = -1;
503 
504  // Keep the points in the same place and move the edge
505  if (hitSurfI == surfI)
506  {
507  newPointi = pointi;
508  }
509  else
510  {
511  newPoints[hitSurfI].append(newPoints[surfI][pointi]);
512  newPointi = newPoints[hitSurfI].size() - 1;
513  }
514 
515  // Split the other face.
516  greenRefine
517  (
518  hitSurf,
519  facei,
520  eIndex,
521  newPointi,
522  newFacesFromSplit
523  );
524 
525  visitedFace[hitSurfI][facei] = true;
526 
527  forAll(newFacesFromSplit, newFacei)
528  {
529  const labelledTri& fN = newFacesFromSplit[newFacei];
530 
531  if (newFacei == 0)
532  {
533  newFaces[hitSurfI][facei] = fN;
534  }
535  else
536  {
537  newFaces[hitSurfI].append(fN);
538  }
539  }
540  }
541  }
542  }
543 
544  Info<< " Number of edges split = " << nChanged << endl;
545 
546  forAll(newSurfaces, surfI)
547  {
548  newSurfaces.set
549  (
550  surfI,
551  new triSurfaceMesh
552  (
553  IOobject
554  (
555  "hookedSurface_" + surfs.names()[surfI],
556  runTime.constant(),
558  runTime
559  ),
560  triSurface
561  (
562  newFaces[surfI],
563  newSurfaces[surfI].patches(),
564  pointField(newPoints[surfI])
565  )
566  )
567  );
568  }
569 
570  } while (nChanged > 0 && nIters <= maxIters);
571 
572  Info<< endl;
573 
574  forAll(newSurfaces, surfI)
575  {
576  const triSurfaceMesh& newSurf = newSurfaces[surfI];
577 
578  Info<< "Writing hooked surface " << newSurf.searchableSurface::name()
579  << endl;
580 
581  newSurf.searchableSurface::write();
582  }
583 
584  Info<< "\nEnd\n" << endl;
585 
586  return 0;
587 }
588 
589 
590 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A bit-packed bool list.
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
bool hit() const
Is there a hit.
Definition: PointHit.H:120
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:54
const Point & hitPoint() const
Return hit point.
const Point & rawPoint() const
Return point with no checking.
label index() const
Return index.
bool hit() const
Is there a hit.
label nEdges() const
Return number of edges in patch.
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const labelList & boundaryPoints() const
Return list of boundary points,.
const labelListList & edgeFaces() const
Return edge-face addressing.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Random number generator.
Definition: Random.H:58
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:74
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:183
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:72
const Type & shapes() const
Reference to shape.
Triangle with additional region number.
Definition: labelledTri.H:60
static const word & geometryDir()
Return the geometry directory name.
Container for searchableSurfaces.
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
treeBoundBox extend(const scalar s) const
Return asymetrically extended bounding box, with guaranteed.
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:54
const labelList & edgeLabels() const
Definition: treeDataEdge.H:179
const edgeList & edges() const
Definition: treeDataEdge.H:169
const pointField & points() const
Definition: treeDataEdge.H:174
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals,...
Triangulated surface description with patch information.
Definition: triSurface.H:69
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label newPointi
Definition: readKivaGrid.H:495
labelList f(nPoints)
dictionary dict
Foam::argList args(argc, argv)
Unit conversion functions.
Random rndGen(label(0))