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 | Copyright (C) 2014-2016 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 
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 //scalar checkEdgeAngle
117 //(
118 // const triSurface& surf,
119 // const label edgeIndex,
120 // const label pointIndex,
121 // const scalar& angle
122 //)
123 //{
124 // const edge& e = surf.edges()[edgeIndex];
125 
126 // vector eVec = e.vec(surf.localPoints());
127 // eVec /= mag(eVec) + SMALL;
128 
129 // const labelList& pEdges = surf.pointEdges()[pointIndex];
130 //
131 // forAll(pEdges, eI)
132 // {
133 // const edge& nearE = surf.edges()[pEdges[eI]];
134 
135 // vector nearEVec = nearE.vec(surf.localPoints());
136 // nearEVec /= mag(nearEVec) + SMALL;
137 
138 // const scalar dot = eVec & nearEVec;
139 // const scalar minCos = degToRad(angle);
140 
141 // if (mag(dot) > minCos)
142 // {
143 // return false;
144 // }
145 // }
146 
147 // return true;
148 //}
149 
150 
151 void createBoundaryEdgeTrees
152 (
153  const PtrList<triSurfaceMesh>& surfs,
155  labelListList& treeBoundaryEdges
156 )
157 {
158  forAll(surfs, surfI)
159  {
160  const triSurface& surf = surfs[surfI];
161 
162  // Boundary edges
163  treeBoundaryEdges[surfI] =
164  labelList
165  (
166  identity(surf.nEdges() - surf.nInternalEdges())
167  + surf.nInternalEdges()
168  );
169 
170  Random rndGen(17301893);
171 
172  // Slightly extended bb. Slightly off-centred just so on symmetric
173  // geometry there are less face/edge aligned items.
174  treeBoundBox bb
175  (
177  );
178 
179  bb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
180  bb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
181 
182  bEdgeTrees.set
183  (
184  surfI,
186  (
188  (
189  false, // cachebb
190  surf.edges(), // edges
191  surf.localPoints(), // points
192  treeBoundaryEdges[surfI] // selected edges
193  ),
194  bb, // bb
195  8, // maxLevel
196  10, // leafsize
197  3.0 // duplicity
198  )
199  );
200  }
201 }
202 
203 
204 class findNearestOpSubset
205 {
206  const indexedOctree<treeDataEdge>& tree_;
207 
208  DynamicList<label>& shapeMask_;
209 
210 public:
211 
212  findNearestOpSubset
213  (
214  const indexedOctree<treeDataEdge>& tree,
215  DynamicList<label>& shapeMask
216  )
217  :
218  tree_(tree),
219  shapeMask_(shapeMask)
220  {}
221 
222  void operator()
223  (
224  const labelUList& indices,
225  const point& sample,
226 
227  scalar& nearestDistSqr,
228  label& minIndex,
229  point& nearestPoint
230  ) const
231  {
232  const treeDataEdge& shape = tree_.shapes();
233 
234  forAll(indices, i)
235  {
236  const label index = indices[i];
237  const label edgeIndex = shape.edgeLabels()[index];
238 
239  if
240  (
241  !shapeMask_.empty()
242  && findIndex(shapeMask_, edgeIndex) != -1
243  )
244  {
245  continue;
246  }
247 
248  const edge& e = shape.edges()[edgeIndex];
249 
250  pointHit nearHit = e.line(shape.points()).nearestDist(sample);
251 
252  // Only register hit if closest point is not an edge point
253  if (nearHit.hit())
254  {
255  scalar distSqr = sqr(nearHit.distance());
256 
257  if (distSqr < nearestDistSqr)
258  {
259  nearestDistSqr = distSqr;
260  minIndex = index;
261  nearestPoint = nearHit.rawPoint();
262  }
263  }
264  }
265  }
266 };
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 int main(int argc, char *argv[])
272 {
274  (
275  "hook surfaces to other surfaces by moving and retriangulating their"
276  "boundary edges to match other surface boundary edges"
277  );
279  argList::validArgs.append("hookTolerance");
280 
281  #include "addDictOption.H"
282 
283  #include "setRootCase.H"
284  #include "createTime.H"
285 
286  const word dictName("surfaceHookUpDict");
288 
289  Info<< "Reading " << dictName << nl << endl;
290 
291  const IOdictionary dict(dictIO);
292 
293  const scalar dist(args.argRead<scalar>(1));
294  const scalar matchTolerance(max(1e-6*dist, SMALL));
295  const label maxIters = 100;
296 
297  Info<< "Hooking distance = " << dist << endl;
298 
299  searchableSurfaces surfs
300  (
301  IOobject
302  (
303  "surfacesToHook",
304  runTime.constant(),
305  "triSurface",
306  runTime
307  ),
308  dict,
309  true // assume single-region names get surface name
310  );
311 
312  Info<< nl << "Reading surfaces: " << endl;
313  forAll(surfs, surfI)
314  {
315  Info<< incrIndent;
316  Info<< nl << indent << "Surface = " << surfs.names()[surfI] << endl;
317 
318  const wordList& regions = surfs[surfI].regions();
319  forAll(regions, surfRegionI)
320  {
321  Info<< incrIndent;
322  Info<< indent << "Regions = " << regions[surfRegionI] << endl;
323  Info<< decrIndent;
324  }
325  Info<< decrIndent;
326  }
327 
328  PtrList<indexedOctree<treeDataEdge>> bEdgeTrees(surfs.size());
329  labelListList treeBoundaryEdges(surfs.size());
330 
331  List<DynamicList<labelledTri>> newFaces(surfs.size());
332  List<DynamicList<point>> newPoints(surfs.size());
333  List<PackedBoolList> visitedFace(surfs.size());
334 
335  PtrList<triSurfaceMesh> newSurfaces(surfs.size());
336  forAll(surfs, surfI)
337  {
338  const triSurfaceMesh& surf =
339  refCast<const triSurfaceMesh>(surfs[surfI]);
340 
341  newSurfaces.set
342  (
343  surfI,
344  new triSurfaceMesh
345  (
346  IOobject
347  (
348  "hookedSurface_" + surfs.names()[surfI],
349  runTime.constant(),
350  "triSurface",
351  runTime
352  ),
353  surf
354  )
355  );
356  }
357 
358  label nChanged = 0;
359  label nIters = 1;
360 
361  do
362  {
363  Info<< nl << "Iteration = " << nIters++ << endl;
364  nChanged = 0;
365 
366  createBoundaryEdgeTrees(newSurfaces, bEdgeTrees, treeBoundaryEdges);
367 
368  forAll(newSurfaces, surfI)
369  {
370  const triSurface& newSurf = newSurfaces[surfI];
371 
372  newFaces[surfI] = newSurf.localFaces();
373  newPoints[surfI] = newSurf.localPoints();
374  visitedFace[surfI] = PackedBoolList(newSurf.size(), false);
375  }
376 
377  forAll(newSurfaces, surfI)
378  {
379  const triSurface& surf = newSurfaces[surfI];
380 
381  List<pointIndexHit> bPointsTobEdges(surf.boundaryPoints().size());
382  labelList bPointsHitTree(surf.boundaryPoints().size(), -1);
383 
384  const labelListList& pointEdges = surf.pointEdges();
385 
386  forAll(bPointsTobEdges, bPointi)
387  {
388  pointIndexHit& nearestHit = bPointsTobEdges[bPointi];
389 
390  const label pointi = surf.boundaryPoints()[bPointi];
391  const point& samplePt = surf.localPoints()[pointi];
392 
393  const labelList& pEdges = pointEdges[pointi];
394 
395  // Add edges connected to the edge to the shapeMask
396  DynamicList<label> shapeMask;
397  shapeMask.append(pEdges);
398 
399  forAll(bEdgeTrees, treeI)
400  {
401  const indexedOctree<treeDataEdge>& bEdgeTree =
402  bEdgeTrees[treeI];
403 
404  pointIndexHit currentHit =
405  bEdgeTree.findNearest
406  (
407  samplePt,
408  sqr(dist),
409  findNearestOpSubset
410  (
411  bEdgeTree,
412  shapeMask
413  )
414  );
415 
416  if
417  (
418  currentHit.hit()
419  &&
420  (
421  !nearestHit.hit()
422  ||
423  (
424  magSqr(currentHit.hitPoint() - samplePt)
425  < magSqr(nearestHit.hitPoint() - samplePt)
426  )
427  )
428  )
429  {
430  nearestHit = currentHit;
431  bPointsHitTree[bPointi] = treeI;
432  }
433  }
434 
435  scalar dist2 = magSqr(nearestHit.rawPoint() - samplePt);
436 
437  if (nearestHit.hit())
438  {
439  // bool rejectEdge =
440  // checkEdgeAngle
441  // (
442  // surf,
443  // nearestHit.index(),
444  // pointi,
445  // 30
446  // );
447 
448  if (dist2 > Foam::sqr(dist))
449  {
450  nearestHit.setMiss();
451  }
452  }
453  }
454 
455  forAll(bPointsTobEdges, bPointi)
456  {
457  const pointIndexHit& eHit = bPointsTobEdges[bPointi];
458 
459  if (eHit.hit())
460  {
461  const label hitSurfI = bPointsHitTree[bPointi];
462  const triSurface& hitSurf = newSurfaces[hitSurfI];
463 
464  const label eIndex =
465  treeBoundaryEdges[hitSurfI][eHit.index()];
466  const edge& e = hitSurf.edges()[eIndex];
467 
468  const label pointi = surf.boundaryPoints()[bPointi];
469 
470  const labelList& eFaces = hitSurf.edgeFaces()[eIndex];
471 
472  if (eFaces.size() != 1)
473  {
475  << "Edge is attached to " << eFaces.size()
476  << " faces." << endl;
477 
478  continue;
479  }
480 
481  const label facei = eFaces[0];
482 
483  if (visitedFace[hitSurfI][facei])
484  {
485  continue;
486  }
487 
488  DynamicList<labelledTri> newFacesFromSplit(2);
489 
490  const point& pt = surf.localPoints()[pointi];
491 
492  if
493  (
494  (
495  magSqr(pt - hitSurf.localPoints()[e.start()])
496  < matchTolerance
497  )
498  || (
499  magSqr(pt - hitSurf.localPoints()[e.end()])
500  < matchTolerance
501  )
502  )
503  {
504  continue;
505  }
506 
507  nChanged++;
508 
509  label newPointi = -1;
510 
511  // Keep the points in the same place and move the edge
512  if (hitSurfI == surfI)
513  {
514  newPointi = pointi;
515  }
516  else
517  {
518  newPoints[hitSurfI].append(newPoints[surfI][pointi]);
519  newPointi = newPoints[hitSurfI].size() - 1;
520  }
521 
522  // Split the other face.
523  greenRefine
524  (
525  hitSurf,
526  facei,
527  eIndex,
528  newPointi,
529  newFacesFromSplit
530  );
531 
532  visitedFace[hitSurfI][facei] = true;
533 
534  forAll(newFacesFromSplit, newFacei)
535  {
536  const labelledTri& fN = newFacesFromSplit[newFacei];
537 
538  if (newFacei == 0)
539  {
540  newFaces[hitSurfI][facei] = fN;
541  }
542  else
543  {
544  newFaces[hitSurfI].append(fN);
545  }
546  }
547  }
548  }
549  }
550 
551  Info<< " Number of edges split = " << nChanged << endl;
552 
553  forAll(newSurfaces, surfI)
554  {
555  newSurfaces.set
556  (
557  surfI,
558  new triSurfaceMesh
559  (
560  IOobject
561  (
562  "hookedSurface_" + surfs.names()[surfI],
563  runTime.constant(),
564  "triSurface",
565  runTime
566  ),
567  triSurface
568  (
569  newFaces[surfI],
570  newSurfaces[surfI].patches(),
571  pointField(newPoints[surfI])
572  )
573  )
574  );
575  }
576 
577  } while (nChanged > 0 && nIters <= maxIters);
578 
579  Info<< endl;
580 
581  forAll(newSurfaces, surfI)
582  {
583  const triSurfaceMesh& newSurf = newSurfaces[surfI];
584 
585  Info<< "Writing hooked surface " << newSurf.searchableSurface::name()
586  << endl;
587 
588  newSurf.searchableSurface::write();
589  }
590 
591  Info<< "\nEnd\n" << endl;
592 
593  return 0;
594 }
595 
596 
597 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
label end() const
Return end vertex label.
Definition: edgeI.H:92
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:169
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
label region() const
Return region label.
Definition: labelledTriI.H:68
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
bool hit() const
Is there a hit.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
const labelList & boundaryPoints() const
Return list of boundary points,.
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Holds data for octree to work on an edges subset.
Definition: treeDataEdge.H:53
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
patches[0]
const edgeList & edges() const
Definition: treeDataEdge.H:169
const labelListList & pointEdges() const
Return point-edge addressing.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
IOoject and searching on triSurface.
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:97
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
Container for searchableSurfaces.
List< label > labelList
A List of labels.
Definition: labelList.H:56
Triangle with additional region number.
Definition: labelledTri.H:57
const Point & rawPoint() const
Return point with no checking.
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: FixedListI.H:115
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Simple random number generator.
Definition: Random.H:49
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label start() const
Return start vertex label.
Definition: edgeI.H:81
static const char nl
Definition: Ostream.H:262
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:177
label nEdges() const
Return number of edges in patch.
const labelList & edgeLabels() const
Definition: treeDataEdge.H:179
A bit-packed bool list.
vector point
Point is a vector.
Definition: point.H:41
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
word dictName("noiseDict")
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const Type & shapes() const
Reference to shape.
label nInternalEdges() const
Number of internal edges.
label index() const
Return index.
const Point & hitPoint() const
Return hit point.
const pointField & points() const
Definition: treeDataEdge.H:174
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
messageStream Info
const labelListList & edgeFaces() const
Return edge-face addressing.
bool hit() const
Is there a hit.
Definition: PointHit.H:120
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
Triangulated surface description with patch information.
Definition: triSurface.H:65
Foam::argList args(argc, argv)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29