surfaceHookUp.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) 2014-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  surfaceHookUp
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
34  Find close open edges and stitches the surface along them
35 
36 Usage
37  - surfaceHookUp hookDistance [OPTION]
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "argList.H"
42 #include "Time.H"
43 
44 #include "triSurfaceMesh.H"
45 #include "indexedOctree.H"
46 #include "treeBoundBox.H"
47 #include "bitSet.H"
48 #include "unitConversion.H"
49 #include "searchableSurfaces.H"
50 #include "IOdictionary.H"
51 
52 using namespace Foam;
53 
54 // Split facei along edgeI at position newPointi
55 void greenRefine
56 (
57  const triSurface& surf,
58  const label facei,
59  const label edgeI,
60  const label newPointi,
61  DynamicList<labelledTri>& newFaces
62 )
63 {
64  const labelledTri& f = surf.localFaces()[facei];
65  const edge& e = surf.edges()[edgeI];
66 
67  // Find index of edge in face.
68 
69  label fp0 = f.find(e[0]);
70  label fp1 = f.fcIndex(fp0);
71  label fp2 = f.fcIndex(fp1);
72 
73  if (f[fp1] == e[1])
74  {
75  // Edge oriented like face
76  newFaces.append
77  (
79  (
80  f[fp0],
81  newPointi,
82  f[fp2],
83  f.region()
84  )
85  );
86  newFaces.append
87  (
89  (
90  newPointi,
91  f[fp1],
92  f[fp2],
93  f.region()
94  )
95  );
96  }
97  else
98  {
99  newFaces.append
100  (
102  (
103  f[fp2],
104  newPointi,
105  f[fp1],
106  f.region()
107  )
108  );
109  newFaces.append
110  (
112  (
113  newPointi,
114  f[fp0],
115  f[fp1],
116  f.region()
117  )
118  );
119  }
120 }
121 
122 
123 //scalar checkEdgeAngle
124 //(
125 // const triSurface& surf,
126 // const label edgeIndex,
127 // const label pointIndex,
128 // const scalar& angle
129 //)
130 //{
131 // const edge& e = surf.edges()[edgeIndex];
132 
133 // const vector eVec = e.unitVec(surf.localPoints());
134 
135 // const labelList& pEdges = surf.pointEdges()[pointIndex];
136 //
137 // forAll(pEdges, eI)
138 // {
139 // const edge& nearE = surf.edges()[pEdges[eI]];
140 
141 // const vector nearEVec = nearE.unitVec(surf.localPoints());
142 
143 // const scalar dot = eVec & nearEVec;
144 // const scalar minCos = degToRad(angle);
145 
146 // if (mag(dot) > minCos)
147 // {
148 // return false;
149 // }
150 // }
151 
152 // return true;
153 //}
154 
155 
156 void createBoundaryEdgeTrees
157 (
158  const PtrList<triSurfaceMesh>& surfs,
160  labelListList& treeBoundaryEdges
161 )
162 {
163  forAll(surfs, surfI)
164  {
165  const triSurface& surf = surfs[surfI];
166 
167  // Boundary edges
168  treeBoundaryEdges[surfI] =
169  identity(surf.nBoundaryEdges(), surf.nInternalEdges());
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  (
177  treeBoundBox(surf.localPoints()).extend(rndGen, 1e-4, ROOTVSMALL)
178  );
179 
180  bEdgeTrees.set
181  (
182  surfI,
184  (
186  (
187  surf.edges(),
188  surf.localPoints(),
189  treeBoundaryEdges[surfI] // selected edges
190  ),
191  bb, // bb
192  8, // maxLevel
193  10, // leafsize
194  3.0 // duplicity
195  )
196  );
197  }
198 }
199 
200 
201 class findNearestOpSubset
202 {
203  const indexedOctree<treeDataEdge>& tree_;
204 
205  DynamicList<label>& shapeMask_;
206 
207 public:
208 
209  findNearestOpSubset
210  (
212  DynamicList<label>& shapeMask
213  )
214  :
215  tree_(tree),
216  shapeMask_(shapeMask)
217  {}
218 
219  void operator()
220  (
221  const labelUList& indices,
222  const point& sample,
223 
224  scalar& nearestDistSqr,
225  label& minIndex,
226  point& nearestPoint
227  ) const
228  {
229  const treeDataEdge& shape = tree_.shapes();
230 
231  for (const label index : indices)
232  {
233  const label edgeIndex = shape.objectIndex(index);
234 
235  if (shapeMask_.found(edgeIndex))
236  {
237  continue;
238  }
239 
240  pointHit nearHit = shape.line(index).nearestDist(sample);
241 
242  // Only register hit if closest point is not an edge point
243  if (nearHit.hit())
244  {
245  const scalar distSqr = sqr(nearHit.distance());
246 
247  if (distSqr < nearestDistSqr)
248  {
249  nearestDistSqr = distSqr;
250  minIndex = index;
251  nearestPoint = nearHit.point();
252  }
253  }
254  }
255  }
256 };
257 
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 int main(int argc, char *argv[])
262 {
264  (
265  "Hook surfaces to other surfaces by moving and retriangulating their"
266  " boundary edges to match other surface boundary edges"
267  );
269  argList::addArgument("hookTolerance", "The point merge tolerance");
270  argList::addOption("dict", "file", "Alternative surfaceHookUpDict");
271 
272  #include "setRootCase.H"
273  #include "createTime.H"
274 
275  const word dictName("surfaceHookUpDict");
277 
278  Info<< "Reading " << dictIO.name() << nl << endl;
279 
280  const IOdictionary dict(dictIO);
281 
282  const scalar dist(args.get<scalar>(1));
283  const scalar matchTolerance(max(1e-6*dist, SMALL));
284  const label maxIters = 100;
285 
286  Info<< "Hooking distance = " << dist << endl;
287 
288  searchableSurfaces surfs
289  (
290  IOobject
291  (
292  "surfacesToHook",
293  runTime.constant(),
294  "triSurface",
295  runTime
296  ),
297  dict,
298  true // assume single-region names get surface name
299  );
300 
301  Info<< nl << "Reading surfaces: " << endl;
302  forAll(surfs, surfI)
303  {
304  Info<< incrIndent;
305  Info<< nl << indent << "Surface = " << surfs.names()[surfI] << endl;
306 
307  const wordList& regions = surfs[surfI].regions();
308  forAll(regions, surfRegionI)
309  {
310  Info<< incrIndent;
311  Info<< indent << "Regions = " << regions[surfRegionI] << endl;
312  Info<< decrIndent;
313  }
314  Info<< decrIndent;
315  }
316 
317  PtrList<indexedOctree<treeDataEdge>> bEdgeTrees(surfs.size());
318  labelListList treeBoundaryEdges(surfs.size());
319 
320  List<DynamicList<labelledTri>> newFaces(surfs.size());
321  List<DynamicList<point>> newPoints(surfs.size());
322  List<bitSet> visitedFace(surfs.size());
323 
324  PtrList<triSurfaceMesh> newSurfaces(surfs.size());
325  forAll(surfs, surfI)
326  {
327  const triSurfaceMesh& surf =
328  refCast<const triSurfaceMesh>(surfs[surfI]);
329 
330  newSurfaces.set
331  (
332  surfI,
333  new triSurfaceMesh
334  (
335  IOobject
336  (
337  "hookedSurface_" + surfs.names()[surfI],
338  runTime.constant(),
339  "triSurface",
340  runTime
341  ),
342  surf
343  )
344  );
345  }
346 
347  label nChanged = 0;
348  label nIters = 1;
349 
350  do
351  {
352  Info<< nl << "Iteration = " << nIters++ << endl;
353  nChanged = 0;
354 
355  createBoundaryEdgeTrees(newSurfaces, bEdgeTrees, treeBoundaryEdges);
356 
357  forAll(newSurfaces, surfI)
358  {
359  const triSurface& newSurf = newSurfaces[surfI];
360 
361  newFaces[surfI] = newSurf.localFaces();
362  newPoints[surfI] = newSurf.localPoints();
363  visitedFace[surfI] = bitSet(newSurf.size(), false);
364  }
365 
366  forAll(newSurfaces, surfI)
367  {
368  const triSurface& surf = newSurfaces[surfI];
369 
370  List<pointIndexHit> bPointsTobEdges(surf.boundaryPoints().size());
371  labelList bPointsHitTree(surf.boundaryPoints().size(), -1);
372 
373  const labelListList& pointEdges = surf.pointEdges();
374 
375  forAll(bPointsTobEdges, bPointi)
376  {
377  pointIndexHit& nearestHit = bPointsTobEdges[bPointi];
378 
379  const label pointi = surf.boundaryPoints()[bPointi];
380  const point& samplePt = surf.localPoints()[pointi];
381 
382  const labelList& pEdges = pointEdges[pointi];
383 
384  // Add edges connected to the edge to the shapeMask
385  DynamicList<label> shapeMask;
386  shapeMask.append(pEdges);
387 
388  forAll(bEdgeTrees, treeI)
389  {
390  const indexedOctree<treeDataEdge>& bEdgeTree =
391  bEdgeTrees[treeI];
392 
393  pointIndexHit currentHit =
394  bEdgeTree.findNearest
395  (
396  samplePt,
397  sqr(dist),
398  findNearestOpSubset
399  (
400  bEdgeTree,
401  shapeMask
402  )
403  );
404 
405  if
406  (
407  currentHit.hit()
408  &&
409  (
410  !nearestHit.hit()
411  ||
412  (
413  currentHit.point().distSqr(samplePt)
414  < nearestHit.point().distSqr(samplePt)
415  )
416  )
417  )
418  {
419  nearestHit = currentHit;
420  bPointsHitTree[bPointi] = treeI;
421  }
422  }
423 
424  if (nearestHit.hit())
425  {
426  // bool rejectEdge =
427  // checkEdgeAngle
428  // (
429  // surf,
430  // nearestHit.index(),
431  // pointi,
432  // 30
433  // );
434 
435  scalar distSqr = nearestHit.point().distSqr(samplePt);
436 
437  if (distSqr > Foam::sqr(dist))
438  {
439  nearestHit.setMiss();
440  }
441  }
442  }
443 
444  forAll(bPointsTobEdges, bPointi)
445  {
446  const pointIndexHit& eHit = bPointsTobEdges[bPointi];
447 
448  if (eHit.hit())
449  {
450  const label hitSurfI = bPointsHitTree[bPointi];
451  const triSurface& hitSurf = newSurfaces[hitSurfI];
452 
453  const label eIndex =
454  treeBoundaryEdges[hitSurfI][eHit.index()];
455  const edge& e = hitSurf.edges()[eIndex];
456 
457  const label pointi = surf.boundaryPoints()[bPointi];
458 
459  const labelList& eFaces = hitSurf.edgeFaces()[eIndex];
460 
461  if (eFaces.size() != 1)
462  {
464  << "Edge is attached to " << eFaces.size()
465  << " faces." << endl;
466 
467  continue;
468  }
469 
470  const label facei = eFaces[0];
471 
472  if (visitedFace[hitSurfI][facei])
473  {
474  continue;
475  }
476 
477  DynamicList<labelledTri> newFacesFromSplit(2);
478 
479  const point& pt = surf.localPoints()[pointi];
480 
481  if
482  (
483  (
484  pt.distSqr(hitSurf.localPoints()[e.start()])
485  < matchTolerance
486  )
487  || (
488  pt.distSqr(hitSurf.localPoints()[e.end()])
489  < matchTolerance
490  )
491  )
492  {
493  continue;
494  }
495 
496  nChanged++;
497 
498  label newPointi = -1;
499 
500  // Keep the points in the same place and move the edge
501  if (hitSurfI == surfI)
502  {
503  newPointi = pointi;
504  }
505  else
506  {
507  newPoints[hitSurfI].append(newPoints[surfI][pointi]);
508  newPointi = newPoints[hitSurfI].size() - 1;
509  }
510 
511  // Split the other face.
512  greenRefine
513  (
514  hitSurf,
515  facei,
516  eIndex,
517  newPointi,
518  newFacesFromSplit
519  );
520 
521  visitedFace[hitSurfI].set(facei);
522 
523  forAll(newFacesFromSplit, newFacei)
524  {
525  const labelledTri& fN = newFacesFromSplit[newFacei];
526 
527  if (newFacei == 0)
528  {
529  newFaces[hitSurfI][facei] = fN;
530  }
531  else
532  {
533  newFaces[hitSurfI].append(fN);
534  }
535  }
536  }
537  }
538  }
539 
540  Info<< " Number of edges split = " << nChanged << endl;
541 
542  forAll(newSurfaces, surfI)
543  {
544  newSurfaces.set
545  (
546  surfI,
547  new triSurfaceMesh
548  (
549  IOobject
550  (
551  "hookedSurface_" + surfs.names()[surfI],
552  runTime.constant(),
553  "triSurface",
554  runTime
555  ),
556  triSurface
557  (
558  newFaces[surfI],
559  newSurfaces[surfI].patches(),
560  pointField(newPoints[surfI])
561  )
562  )
563  );
564  }
565 
566  } while (nChanged > 0 && nIters <= maxIters);
567 
568  Info<< endl;
569 
570  forAll(newSurfaces, surfI)
571  {
572  const triSurfaceMesh& newSurf = newSurfaces[surfI];
573 
574  Info<< "Writing hooked surface " << newSurf.searchableSurface::name()
575  << endl;
576 
577  newSurf.searchableSurface::write();
578  }
579 
580  Info<< "\nEnd\n" << endl;
581 
582  return 0;
583 }
584 
585 
586 // ************************************************************************* //
const labelListList & pointEdges() const
Return point-edge addressing.
void setMiss() noexcept
Set the hit status off.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
dictionary dict
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Type & shapes() const noexcept
Reference to shape.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
const word dictName("faMeshDefinition")
const labelList & boundaryPoints() const
Return list of boundary points, address into LOCAL point list.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool found(const T &val, label pos=0) const
Same as contains()
Definition: UList.H:879
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
label nInternalEdges() const
Number of internal edges.
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
label objectIndex(const label index) const
Map from shape index to original (non-subset) edge label.
Definition: treeDataEdge.H:388
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
IOoject and searching on triSurface.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
bool hit() const noexcept
Is there a hit.
Definition: pointHit.H:145
const point_type & point() const noexcept
Return point, no checks.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
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
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Tree tree(triangles.begin(), triangles.end())
Container for searchableSurfaces. The collection is specified as a dictionary. For example...
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
A triFace with additional (region) index.
Definition: labelledTri.H:53
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
label nBoundaryEdges() const
Number of boundary edges == (nEdges() - nInternalEdges())
label index() const noexcept
Return the hit index.
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
Random number generator.
Definition: Random.H:55
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
labelList f(nPoints)
bool hit() const noexcept
Is there a hit?
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
Non-pointer based hierarchical recursive searching.
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.
label newPointi
Definition: readKivaGrid.H:496
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
scalar distance() const noexcept
Return distance to hit.
Definition: pointHit.H:169
const polyBoundaryMesh & patches
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
const point_type & point() const noexcept
Return the point, no checks.
Definition: pointHit.H:161
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
static constexpr int maxIters
const linePointRef line(const label index) const
Geometric line for edge at specified shape index. Frequently used.
Definition: treeDataEdge.H:404
Triangulated surface description with patch information.
Definition: triSurface.H:71
Foam::argList args(argc, argv)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:167
Namespace for OpenFOAM.