Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 surfaceHookUp
29
30Group
31 grpSurfaceUtilities
32
33Description
34 Find close open edges and stitches the surface along them
35
36Usage
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
52using namespace Foam;
53
54// Split facei along edgeI at position newPointi
55void greenRefine
56(
57 const triSurface& surf,
58 const label facei,
59 const label edgeI,
60 const label newPointi,
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],
82 f[fp2],
83 f.region()
84 )
85 );
86 newFaces.append
87 (
89 (
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
156void 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
201class findNearestOpSubset
202{
203 const indexedOctree<treeDataEdge>& tree_;
204
205 DynamicList<label>& shapeMask_;
206
207public:
208
209 findNearestOpSubset
210 (
211 const indexedOctree<treeDataEdge>& tree,
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
261int 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");
272 (
273 "maxIters",
274 "number",
275 "Maximum number of iterations (default: 100)"
276 );
277
278 #include "setRootCase.H"
279 #include "createTime.H"
280
281 const word dictName("surfaceHookUpDict");
283
284 Info<< "Reading " << dictIO.name() << nl << endl;
285
286 const IOdictionary dict(dictIO);
287
288 const scalar dist(args.get<scalar>(1));
289 const scalar matchTolerance(Foam::max(1e-6*dist, SMALL));
290 const label maxIters = args.getOrDefault<label>("maxIters", 100);
291
292 Info<< "Hooking distance = " << dist << nl
293 << "Maximum iterations = " << maxIters << endl;
294
296 (
298 (
299 "surfacesToHook",
300 runTime.constant(),
301 "triSurface",
302 runTime
303 ),
304 dict,
305 true // assume single-region names get surface name
306 );
307
308 Info<< nl << "Reading surfaces: " << endl;
309 forAll(surfs, surfI)
310 {
312 Info<< nl << indent << "Surface = " << surfs.names()[surfI] << endl;
313
314 const wordList& regions = surfs[surfI].regions();
315 forAll(regions, surfRegionI)
316 {
318 Info<< indent << "Regions = " << regions[surfRegionI] << endl;
320 }
322 }
323
324 PtrList<indexedOctree<treeDataEdge>> bEdgeTrees(surfs.size());
325 labelListList treeBoundaryEdges(surfs.size());
326
327 List<DynamicList<labelledTri>> newFaces(surfs.size());
328 List<DynamicList<point>> newPoints(surfs.size());
329 List<bitSet> visitedFace(surfs.size());
330
331 PtrList<triSurfaceMesh> newSurfaces(surfs.size());
332 forAll(surfs, surfI)
333 {
334 const triSurfaceMesh& surf =
335 refCast<const triSurfaceMesh>(surfs[surfI]);
336
337 newSurfaces.set
338 (
339 surfI,
341 (
343 (
344 "hookedSurface_" + surfs.names()[surfI],
345 runTime.constant(),
346 "triSurface",
347 runTime
348 ),
349 surf
350 )
351 );
352 }
353
354 label nChanged = 0;
355 label nIters = 1;
356
357 do
358 {
359 Info<< nl << "Iteration = " << nIters++ << endl;
360 nChanged = 0;
361
362 createBoundaryEdgeTrees(newSurfaces, bEdgeTrees, treeBoundaryEdges);
363
364 forAll(newSurfaces, surfI)
365 {
366 const triSurface& newSurf = newSurfaces[surfI];
367
368 newFaces[surfI] = newSurf.localFaces();
369 newPoints[surfI] = newSurf.localPoints();
370 visitedFace[surfI] = bitSet(newSurf.size(), false);
371 }
372
373 forAll(newSurfaces, surfI)
374 {
375 const triSurface& surf = newSurfaces[surfI];
376
377 List<pointIndexHit> bPointsTobEdges(surf.boundaryPoints().size());
378 labelList bPointsHitTree(surf.boundaryPoints().size(), -1);
379
380 const labelListList& pointEdges = surf.pointEdges();
381
382 forAll(bPointsTobEdges, bPointi)
383 {
384 pointIndexHit& nearestHit = bPointsTobEdges[bPointi];
385
386 const label pointi = surf.boundaryPoints()[bPointi];
387 const point& samplePt = surf.localPoints()[pointi];
388
389 const labelList& pEdges = pointEdges[pointi];
390
391 // Add edges connected to the edge to the shapeMask
392 DynamicList<label> shapeMask;
393 shapeMask.append(pEdges);
394
395 forAll(bEdgeTrees, treeI)
396 {
397 const indexedOctree<treeDataEdge>& bEdgeTree =
398 bEdgeTrees[treeI];
399
400 pointIndexHit currentHit =
401 bEdgeTree.findNearest
402 (
403 samplePt,
404 sqr(dist),
405 findNearestOpSubset
406 (
407 bEdgeTree,
408 shapeMask
409 )
410 );
411
412 if
413 (
414 currentHit.hit()
415 &&
416 (
417 !nearestHit.hit()
418 ||
419 (
420 currentHit.point().distSqr(samplePt)
421 < nearestHit.point().distSqr(samplePt)
422 )
423 )
424 )
425 {
426 nearestHit = currentHit;
427 bPointsHitTree[bPointi] = treeI;
428 }
429 }
430
431 if (nearestHit.hit())
432 {
433 // bool rejectEdge =
434 // checkEdgeAngle
435 // (
436 // surf,
437 // nearestHit.index(),
438 // pointi,
439 // 30
440 // );
441
442 scalar distSqr = nearestHit.point().distSqr(samplePt);
443
444 if (distSqr > Foam::sqr(dist))
445 {
446 nearestHit.setMiss();
447 }
448 }
449 }
450
451 forAll(bPointsTobEdges, bPointi)
452 {
453 const pointIndexHit& eHit = bPointsTobEdges[bPointi];
454
455 if (eHit.hit())
456 {
457 const label hitSurfI = bPointsHitTree[bPointi];
458 const triSurface& hitSurf = newSurfaces[hitSurfI];
459
460 const label eIndex =
461 treeBoundaryEdges[hitSurfI][eHit.index()];
462 const edge& e = hitSurf.edges()[eIndex];
463
464 const label pointi = surf.boundaryPoints()[bPointi];
465
466 const labelList& eFaces = hitSurf.edgeFaces()[eIndex];
467
468 if (eFaces.size() != 1)
469 {
471 << "Edge is attached to " << eFaces.size()
472 << " faces." << endl;
473
474 continue;
475 }
476
477 const label facei = eFaces[0];
478
479 if (visitedFace[hitSurfI][facei])
480 {
481 continue;
482 }
483
484 DynamicList<labelledTri> newFacesFromSplit(2);
485
486 const point& pt = surf.localPoints()[pointi];
487
488 if
489 (
490 (
491 pt.distSqr(hitSurf.localPoints()[e.start()])
492 < matchTolerance
493 )
494 || (
495 pt.distSqr(hitSurf.localPoints()[e.end()])
496 < matchTolerance
497 )
498 )
499 {
500 continue;
501 }
502
503 nChanged++;
504
505 label newPointi = -1;
506
507 // Keep the points in the same place and move the edge
508 if (hitSurfI == surfI)
509 {
510 newPointi = pointi;
511 }
512 else
513 {
514 newPoints[hitSurfI].append(newPoints[surfI][pointi]);
515 newPointi = newPoints[hitSurfI].size() - 1;
516 }
517
518 // Split the other face.
519 greenRefine
520 (
521 hitSurf,
522 facei,
523 eIndex,
524 newPointi,
525 newFacesFromSplit
526 );
527
528 visitedFace[hitSurfI].set(facei);
529
530 forAll(newFacesFromSplit, newFacei)
531 {
532 const labelledTri& fN = newFacesFromSplit[newFacei];
533
534 if (newFacei == 0)
535 {
536 newFaces[hitSurfI][facei] = fN;
537 }
538 else
539 {
540 newFaces[hitSurfI].append(fN);
541 }
542 }
543 }
544 }
545 }
546
547 Info<< " Number of edges split = " << nChanged << endl;
548
549 forAll(newSurfaces, surfI)
550 {
551 newSurfaces.set
552 (
553 surfI,
555 (
557 (
558 "hookedSurface_" + surfs.names()[surfI],
559 runTime.constant(),
560 "triSurface",
561 runTime
562 ),
564 (
565 newFaces[surfI],
566 newSurfaces[surfI].patches(),
567 pointField(newPoints[surfI])
568 )
569 )
570 );
571 }
572
573 } while (nChanged > 0 && nIters <= maxIters);
574
575 Info<< endl;
576
577 forAll(newSurfaces, surfI)
578 {
579 const triSurfaceMesh& newSurf = newSurfaces[surfI];
580
581 Info<< "Writing hooked surface " << newSurf.searchableSurface::name()
582 << endl;
583
584 newSurf.searchableSurface::write();
585 }
586
587 Info<< "\nEnd\n" << endl;
588
589 return 0;
590}
591
592
593// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
bool hit() const noexcept
Is there a hit.
Definition pointHit.H:145
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
label index() const noexcept
Return the hit index.
void setMiss() noexcept
Set the hit status off.
bool hit() const noexcept
Is there a hit?
const point_type & point() const noexcept
Return point, no checks.
label nBoundaryEdges() const
Number of boundary edges == (nEdges() - nInternalEdges()).
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 Field< point_type > & localPoints() const
Return pointField of points in patch.
const labelList & boundaryPoints() const
Return list of boundary points, address into LOCAL point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Random number generator.
Definition Random.H:56
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition argList.C:366
static void noParallel()
Remove the parallel options.
Definition argList.C:599
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:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
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:62
Non-pointer based hierarchical recursive searching.
const Type & shapes() const noexcept
Reference to shape.
A triFace with additional (region) index.
Definition labelledTri.H:56
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition lineI.H:180
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Standard boundBox with extra functionality for use in octree.
treeBoundBox extend(Random &rndGen, const scalar s) const
Return slightly wider bounding box.
Holds data for octree to work on an edges subset.
const linePointRef line(const label index) const
Geometric line for edge at specified shape index. Frequently used.
label objectIndex(const label index) const
Map from shape index to original (non-subset) edge label.
IOoject and searching on triSurface.
Triangulated surface description with patch information.
Definition triSurface.H:74
A class for handling words, derived from Foam::string.
Definition word.H:66
const polyBoundaryMesh & patches
engineTime & runTime
const word dictName("faMeshDefinition")
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
static constexpr int maxIters
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...
vector point
Point is a vector.
Definition point.H:37
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label newPointi
labelList f(nPoints)
dictionary dict
Tree tree(triangles.begin(), triangles.end())
IOobject dictIO
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.
Random rndGen