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-2024 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 "searchableSurfaces.H"
43 #include "systemDict.H"
44 
45 using namespace Foam;
46 
47 // Split facei along edgeI at position newPointi
48 void greenRefine
49 (
50  const triSurface& surf,
51  const label facei,
52  const label edgeI,
53  const label newPointi,
54  DynamicList<labelledTri>& newFaces
55 )
56 {
57  const labelledTri& f = surf.localFaces()[facei];
58  const edge& e = surf.edges()[edgeI];
59 
60  // Find index of edge in face.
61 
62  label fp0 = findIndex(f, e[0]);
63  label fp1 = f.fcIndex(fp0);
64  label fp2 = f.fcIndex(fp1);
65 
66  if (f[fp1] == e[1])
67  {
68  // Edge oriented like face
69  newFaces.append
70  (
72  (
73  f[fp0],
74  newPointi,
75  f[fp2],
76  f.region()
77  )
78  );
79  newFaces.append
80  (
82  (
83  newPointi,
84  f[fp1],
85  f[fp2],
86  f.region()
87  )
88  );
89  }
90  else
91  {
92  newFaces.append
93  (
95  (
96  f[fp2],
97  newPointi,
98  f[fp1],
99  f.region()
100  )
101  );
102  newFaces.append
103  (
105  (
106  newPointi,
107  f[fp0],
108  f[fp1],
109  f.region()
110  )
111  );
112  }
113 }
114 
115 
116 void createBoundaryEdgeTrees
117 (
118  const PtrList<triSurfaceMesh>& surfs,
120  labelListList& treeBoundaryEdges
121 )
122 {
123  forAll(surfs, surfI)
124  {
125  const triSurface& surf = surfs[surfI];
126 
127  // Boundary edges
128  treeBoundaryEdges[surfI] =
129  labelList
130  (
131  identityMap(surf.nEdges() - surf.nInternalEdges())
132  + surf.nInternalEdges()
133  );
134 
135  randomGenerator rndGen(17301893);
136 
137  // Slightly extended bb. Slightly off-centred just so on symmetric
138  // geometry there are less face/edge aligned items.
139  treeBoundBox bb
140  (
142  );
143 
144  bEdgeTrees.set
145  (
146  surfI,
148  (
150  (
151  false, // cachebb
152  surf.edges(), // edges
153  surf.localPoints(), // points
154  treeBoundaryEdges[surfI] // selected edges
155  ),
156  bb, // bb
157  8, // maxLevel
158  10, // leafsize
159  3.0 // duplicity
160  )
161  );
162  }
163 }
164 
165 
166 class findNearestOpSubset
167 {
168  const indexedOctree<treeDataEdge>& tree_;
169 
170  DynamicList<label>& shapeMask_;
171 
172 public:
173 
174  findNearestOpSubset
175  (
176  const indexedOctree<treeDataEdge>& tree,
177  DynamicList<label>& shapeMask
178  )
179  :
180  tree_(tree),
181  shapeMask_(shapeMask)
182  {}
183 
184  void operator()
185  (
186  const labelUList& indices,
187  const point& sample,
188 
189  scalar& nearestDistSqr,
190  label& minIndex,
191  point& nearestPoint
192  ) const
193  {
194  const treeDataEdge& shape = tree_.shapes();
195 
196  forAll(indices, i)
197  {
198  const label index = indices[i];
199  const label edgeIndex = shape.edgeLabels()[index];
200 
201  if
202  (
203  !shapeMask_.empty()
204  && findIndex(shapeMask_, edgeIndex) != -1
205  )
206  {
207  continue;
208  }
209 
210  const edge& e = shape.edges()[edgeIndex];
211 
212  pointHit nearHit = e.line(shape.points()).nearestDist(sample);
213 
214  // Only register hit if closest point is not an edge point
215  if (nearHit.hit())
216  {
217  scalar distSqr = sqr(nearHit.distance());
218 
219  if (distSqr < nearestDistSqr)
220  {
221  nearestDistSqr = distSqr;
222  minIndex = index;
223  nearestPoint = nearHit.rawPoint();
224  }
225  }
226  }
227  }
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 int main(int argc, char *argv[])
234 {
236  (
237  "hook surfaces to other surfaces by moving and retriangulating their"
238  "boundary edges to match other surface boundary edges"
239  );
241  argList::validArgs.append("hookTolerance");
242 
243  #include "addDictOption.H"
244 
245  #include "setRootCase.H"
246  #include "createTime.H"
247 
248  const dictionary dict(systemDict("surfaceHookUpDict", args, runTime));
249 
250  const scalar dist(args.argRead<scalar>(1));
251  const scalar matchTolerance(max(1e-6*dist, small));
252  const label maxIters = 100;
253 
254  Info<< "Hooking distance = " << dist << endl;
255 
256  searchableSurfaces surfs
257  (
258  IOobject
259  (
260  "surfacesToHook",
261  runTime.constant(),
263  runTime
264  ),
265  dict,
266  true // assume single-region names get surface name
267  );
268 
269  Info<< nl << "Reading surfaces: " << endl;
270  forAll(surfs, surfI)
271  {
272  Info<< incrIndent;
273  Info<< nl << indent << "Surface = " << surfs.names()[surfI] << endl;
274 
275  const wordList& regions = surfs[surfI].regions();
276  forAll(regions, surfRegionI)
277  {
278  Info<< incrIndent;
279  Info<< indent << "Regions = " << regions[surfRegionI] << endl;
280  Info<< decrIndent;
281  }
282  Info<< decrIndent;
283  }
284 
285  PtrList<indexedOctree<treeDataEdge>> bEdgeTrees(surfs.size());
286  labelListList treeBoundaryEdges(surfs.size());
287 
288  List<DynamicList<labelledTri>> newFaces(surfs.size());
289  List<DynamicList<point>> newPoints(surfs.size());
290  List<PackedBoolList> visitedFace(surfs.size());
291 
292  PtrList<triSurfaceMesh> newSurfaces(surfs.size());
293  forAll(surfs, surfI)
294  {
295  const triSurfaceMesh& surf =
296  refCast<const triSurfaceMesh>(surfs[surfI]);
297 
298  newSurfaces.set
299  (
300  surfI,
301  new triSurfaceMesh
302  (
303  IOobject
304  (
305  "hookedSurface_" + surfs.names()[surfI],
306  runTime.constant(),
308  runTime
309  ),
310  surf
311  )
312  );
313  }
314 
315  label nChanged = 0;
316  label nIters = 1;
317 
318  do
319  {
320  Info<< nl << "Iteration = " << nIters++ << endl;
321  nChanged = 0;
322 
323  createBoundaryEdgeTrees(newSurfaces, bEdgeTrees, treeBoundaryEdges);
324 
325  forAll(newSurfaces, surfI)
326  {
327  const triSurface& newSurf = newSurfaces[surfI];
328 
329  newFaces[surfI] = newSurf.localFaces();
330  newPoints[surfI] = newSurf.localPoints();
331  visitedFace[surfI] = PackedBoolList(newSurf.size(), false);
332  }
333 
334  forAll(newSurfaces, surfI)
335  {
336  const triSurface& surf = newSurfaces[surfI];
337 
338  List<pointIndexHit> bPointsTobEdges(surf.boundaryPoints().size());
339  labelList bPointsHitTree(surf.boundaryPoints().size(), -1);
340 
341  const labelListList& pointEdges = surf.pointEdges();
342 
343  forAll(bPointsTobEdges, bPointi)
344  {
345  pointIndexHit& nearestHit = bPointsTobEdges[bPointi];
346 
347  const label pointi = surf.boundaryPoints()[bPointi];
348  const point& samplePt = surf.localPoints()[pointi];
349 
350  const labelList& pEdges = pointEdges[pointi];
351 
352  // Add edges connected to the edge to the shapeMask
353  DynamicList<label> shapeMask;
354  shapeMask.append(pEdges);
355 
356  forAll(bEdgeTrees, treeI)
357  {
358  const indexedOctree<treeDataEdge>& bEdgeTree =
359  bEdgeTrees[treeI];
360 
361  pointIndexHit currentHit =
362  bEdgeTree.findNearest
363  (
364  samplePt,
365  sqr(dist),
366  findNearestOpSubset
367  (
368  bEdgeTree,
369  shapeMask
370  )
371  );
372 
373  if
374  (
375  currentHit.hit()
376  &&
377  (
378  !nearestHit.hit()
379  ||
380  (
381  magSqr(currentHit.hitPoint() - samplePt)
382  < magSqr(nearestHit.hitPoint() - samplePt)
383  )
384  )
385  )
386  {
387  nearestHit = currentHit;
388  bPointsHitTree[bPointi] = treeI;
389  }
390  }
391 
392  scalar dist2 = magSqr(nearestHit.rawPoint() - samplePt);
393 
394  if (nearestHit.hit())
395  {
396  if (dist2 > Foam::sqr(dist))
397  {
398  nearestHit.setMiss();
399  }
400  }
401  }
402 
403  forAll(bPointsTobEdges, bPointi)
404  {
405  const pointIndexHit& eHit = bPointsTobEdges[bPointi];
406 
407  if (eHit.hit())
408  {
409  const label hitSurfI = bPointsHitTree[bPointi];
410  const triSurface& hitSurf = newSurfaces[hitSurfI];
411 
412  const label eIndex =
413  treeBoundaryEdges[hitSurfI][eHit.index()];
414  const edge& e = hitSurf.edges()[eIndex];
415 
416  const label pointi = surf.boundaryPoints()[bPointi];
417 
418  const labelList& eFaces = hitSurf.edgeFaces()[eIndex];
419 
420  if (eFaces.size() != 1)
421  {
423  << "Edge is attached to " << eFaces.size()
424  << " faces." << endl;
425 
426  continue;
427  }
428 
429  const label facei = eFaces[0];
430 
431  if (visitedFace[hitSurfI][facei])
432  {
433  continue;
434  }
435 
436  DynamicList<labelledTri> newFacesFromSplit(2);
437 
438  const point& pt = surf.localPoints()[pointi];
439 
440  if
441  (
442  (
443  magSqr(pt - hitSurf.localPoints()[e.start()])
444  < matchTolerance
445  )
446  || (
447  magSqr(pt - hitSurf.localPoints()[e.end()])
448  < matchTolerance
449  )
450  )
451  {
452  continue;
453  }
454 
455  nChanged++;
456 
457  label newPointi = -1;
458 
459  // Keep the points in the same place and move the edge
460  if (hitSurfI == surfI)
461  {
462  newPointi = pointi;
463  }
464  else
465  {
466  newPoints[hitSurfI].append(newPoints[surfI][pointi]);
467  newPointi = newPoints[hitSurfI].size() - 1;
468  }
469 
470  // Split the other face.
471  greenRefine
472  (
473  hitSurf,
474  facei,
475  eIndex,
476  newPointi,
477  newFacesFromSplit
478  );
479 
480  visitedFace[hitSurfI][facei] = true;
481 
482  forAll(newFacesFromSplit, newFacei)
483  {
484  const labelledTri& fN = newFacesFromSplit[newFacei];
485 
486  if (newFacei == 0)
487  {
488  newFaces[hitSurfI][facei] = fN;
489  }
490  else
491  {
492  newFaces[hitSurfI].append(fN);
493  }
494  }
495  }
496  }
497  }
498 
499  Info<< " Number of edges split = " << nChanged << endl;
500 
501  forAll(newSurfaces, surfI)
502  {
503  newSurfaces.set
504  (
505  surfI,
506  new triSurfaceMesh
507  (
508  IOobject
509  (
510  "hookedSurface_" + surfs.names()[surfI],
511  runTime.constant(),
513  runTime
514  ),
515  triSurface
516  (
517  newFaces[surfI],
518  newSurfaces[surfI].patches(),
519  pointField(newPoints[surfI])
520  )
521  )
522  );
523  }
524 
525  } while (nChanged > 0 && nIters <= maxIters);
526 
527  Info<< endl;
528 
529  forAll(newSurfaces, surfI)
530  {
531  const triSurfaceMesh& newSurf = newSurfaces[surfI];
532 
533  Info<< "Writing hooked surface " << newSurf.searchableSurface::name()
534  << endl;
535 
536  newSurf.searchableSurface::write();
537  }
538 
539  Info<< "\nEnd\n" << endl;
540 
541  return 0;
542 }
543 
544 
545 // ************************************************************************* //
#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
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:162
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
Random number generator.
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 asymmetrically 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:106
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
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:257
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
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:227
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label newPointi
Definition: readKivaGrid.H:495
labelList f(nPoints)
dictionary dict
randomGenerator rndGen(653213)
Foam::argList args(argc, argv)