slidingInterfaceProjectPoints.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 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 \*---------------------------------------------------------------------------*/
28 
29 #include "slidingInterface.H"
30 #include "polyMesh.H"
31 #include "line.H"
32 #include "polyTopoChanger.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 const Foam::scalar Foam::slidingInterface::pointMergeTolDefault_ = 0.05;
37 const Foam::scalar Foam::slidingInterface::edgeMergeTolDefault_ = 0.01;
38 const Foam::label Foam::slidingInterface::nFacesPerSlaveEdgeDefault_ = 5;
39 const Foam::label Foam::slidingInterface::edgeFaceEscapeLimitDefault_ = 10;
40 
41 const Foam::scalar Foam::slidingInterface::integralAdjTolDefault_ = 0.05;
42 const Foam::scalar
43  Foam::slidingInterface::edgeMasterCatchFractionDefault_ = 0.4;
44 const Foam::scalar Foam::slidingInterface::edgeEndCutoffTolDefault_ = 0.0001;
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 // Index of debug signs:
50 // a - integral match adjustment: point adjusted
51 // n - integral match adjustment: point not adjusted
52 // . - loop of the edge-to-face interaction detection
53 // x - reversal of direction in edge-to-face interaction detection
54 // + - complete edge-to-face interaction detection
55 // z - incomplete edge-to-face interaction detection. This may be OK for edges
56 // crossing from one to the other side of multiply connected master patch
57 // * - colinear triangle: adjusting projection with slave face normal
58 // m - master point inserted into the edge
59 
60 bool Foam::slidingInterface::projectPoints() const
61 {
62  if (debug)
63  {
64  Pout<< FUNCTION_NAME << nl
65  << ": for object " << name() << " : "
66  << "Projecting slave points onto master surface." << endl;
67  }
68 
69  // Algorithm:
70  // 1) Go through all the points of the master and slave patch and calculate
71  // minimum edge length coming from the point. Calculate the point
72  // merge tolerance as the fraction of minimum edge length.
73  // 2) Project all the slave points onto the master patch
74  // in the normal direction.
75  // 3) If some points have missed and the match is integral, the
76  // points in question will be adjusted. Find the nearest face for
77  // those points and continue with the procedure.
78  // 4) For every point, check all the points of the face it has hit.
79  // For every pair of points find if their distance is smaller than
80  // both the master and slave merge tolerance. If so, the slave point
81  // is moved to the location of the master point. Remember the master
82  // point index.
83  // 5) For every unmerged slave point, check its distance to all the
84  // edges of the face it has hit. If the distance is smaller than the
85  // edge merge tolerance, the point will be moved onto the edge.
86  // Remember the master edge index.
87  // 6) The remaning slave points will be projected into faces. Remember the
88  // master face index.
89  // 7) For the points that miss the master patch, grab the nearest face
90  // on the master and leave the slave point where it started
91  // from and the miss is recorded.
92 
93  const polyMesh& mesh = topoChanger().mesh();
94 
95  const primitiveFacePatch& masterPatch =
96  mesh.faceZones()[masterFaceZoneID_.index()]();
97 
98  const primitiveFacePatch& slavePatch =
99  mesh.faceZones()[slaveFaceZoneID_.index()]();
100 
101  // Get references to local points, local edges and local faces
102  // for master and slave patch
103 // const labelList& masterMeshPoints = masterPatch.meshPoints();
104  const pointField& masterLocalPoints = masterPatch.localPoints();
105  const faceList& masterLocalFaces = masterPatch.localFaces();
106  const edgeList& masterEdges = masterPatch.edges();
107  const labelListList& masterEdgeFaces = masterPatch.edgeFaces();
108  const labelListList& masterFaceEdges = masterPatch.faceEdges();
109  const labelListList& masterFaceFaces = masterPatch.faceFaces();
110 // Pout<< "Master patch. Local points: " << masterLocalPoints << nl
111 // << "Master patch. Mesh points: " << masterPatch.meshPoints() << nl
112 // << "Local faces: " << masterLocalFaces << nl
113 // << "local edges: " << masterEdges << endl;
114 
115 // const labelList& slaveMeshPoints = slavePatch.meshPoints();
116  const pointField& slaveLocalPoints = slavePatch.localPoints();
117  const edgeList& slaveEdges = slavePatch.edges();
118  const labelListList& slaveEdgeFaces = slavePatch.edgeFaces();
119  const vectorField& slavePointNormals = slavePatch.pointNormals();
120 // const vectorField& slaveFaceNormals = slavePatch.faceNormals();
121 // Pout<< "Slave patch. Local points: " << slaveLocalPoints << nl
122 // << "Slave patch. Mesh points: " << slavePatch.meshPoints() << nl
123 // << "Local faces: " << slavePatch.localFaces() << nl
124 // << "local edges: " << slaveEdges << endl;
125 
126  // Calculate min edge distance for points and faces
127 
128  // Calculate min edge length for the points and faces of master patch
129  scalarField minMasterPointLength(masterLocalPoints.size(), GREAT);
130  scalarField minMasterFaceLength(masterPatch.size(), GREAT);
131 
132  forAll(masterEdges, edgei)
133  {
134  const edge& curEdge = masterEdges[edgei];
135 
136  const scalar curLength = curEdge.mag(masterLocalPoints);
137 
138  // Do points
139  minMasterPointLength[curEdge.start()] =
140  min
141  (
142  minMasterPointLength[curEdge.start()],
143  curLength
144  );
145 
146  minMasterPointLength[curEdge.end()] =
147  min
148  (
149  minMasterPointLength[curEdge.end()],
150  curLength
151  );
152 
153  // Do faces
154  const labelList& curFaces = masterEdgeFaces[edgei];
155 
156  for (const label facei : curFaces)
157  {
158  minMasterFaceLength[facei] =
159  min
160  (
161  minMasterFaceLength[facei],
162  curLength
163  );
164  }
165  }
166 
167 // Pout<< "min length for master points: " << minMasterPointLength << endl
168 // << "min length for master faces: " << minMasterFaceLength << endl;
169 
170  // Calculate min edge length for the points and faces of slave patch
171  scalarField minSlavePointLength(slaveLocalPoints.size(), GREAT);
172  scalarField minSlaveFaceLength(slavePatch.size(), GREAT);
173 
174  forAll(slaveEdges, edgei)
175  {
176  const edge& curEdge = slaveEdges[edgei];
177 
178  const scalar curLength = curEdge.mag(slaveLocalPoints);
179 
180  // Do points
181  minSlavePointLength[curEdge.start()] =
182  min
183  (
184  minSlavePointLength[curEdge.start()],
185  curLength
186  );
187 
188  minSlavePointLength[curEdge.end()] =
189  min
190  (
191  minSlavePointLength[curEdge.end()],
192  curLength
193  );
194 
195  // Do faces
196  const labelList& curFaces = slaveEdgeFaces[edgei];
197 
198  for (const label facei : curFaces)
199  {
200  minSlaveFaceLength[facei] =
201  min
202  (
203  minSlaveFaceLength[facei],
204  curLength
205  );
206  }
207  }
208 
209 // Pout<< "min length for slave points: " << minSlavePointLength << endl
210 // << "min length for slave faces: " << minSlaveFaceLength << endl;
211 
212  // Project slave points onto the master patch
213 
214  // Face hit by the slave point
215  List<objectHit> slavePointFaceHits =
216  slavePatch.projectPoints
217  (
218  masterPatch,
219  slavePointNormals,
220  projectionAlgo_
221  );
222 
223  if (debug)
224  {
225  label nHits = 0;
226 
227  for (const auto& hit : slavePointFaceHits)
228  {
229  if (hit.hit())
230  {
231  ++nHits;
232  }
233  }
234 
235  Pout<< "Number of hits in point projection: " << nHits
236  << " out of " << slavePointFaceHits.size() << " points."
237  << endl;
238  }
239 
240  // Projected slave points are stored for mesh motion correction
241  projectedSlavePointsPtr_.reset
242  (
243  new pointField(slavePointFaceHits.size(), Zero)
244  );
245  auto& projectedSlavePoints = *projectedSlavePointsPtr_;
246 
247  // Adjust projection to type of match
248 
249  label nAdjustedPoints = 0;
250 
251  // If the match is integral and some points have missed,
252  // find nearest master face
253  if (matchType_ == INTEGRAL)
254  {
255  if (debug)
256  {
257  Pout<< "bool slidingInterface::projectPoints() for object "
258  << name() << " : "
259  << "Adjusting point projection for integral match: ";
260  }
261 
262  forAll(slavePointFaceHits, pointi)
263  {
264  if (slavePointFaceHits[pointi].hit())
265  {
266  // Grab the hit point
267  projectedSlavePoints[pointi] =
268  masterLocalFaces
269  [slavePointFaceHits[pointi].hitObject()].ray
270  (
271  slaveLocalPoints[pointi],
272  slavePointNormals[pointi],
273  masterLocalPoints,
274  projectionAlgo_
275  ).hitPoint();
276  }
277  else
278  {
279  // Grab the nearest point on the face (edge)
280  pointHit missAdjust =
281  masterLocalFaces[slavePointFaceHits[pointi].hitObject()].ray
282  (
283  slaveLocalPoints[pointi],
284  slavePointNormals[pointi],
285  masterLocalPoints,
286  projectionAlgo_
287  );
288 
289  const point nearPoint = missAdjust.missPoint();
290  const point missPoint =
291  slaveLocalPoints[pointi]
292  + slavePointNormals[pointi]*missAdjust.distance();
293 
294  // Calculate the tolerance
295  const scalar mergeTol =
296  integralAdjTol_*minSlavePointLength[pointi];
297 
298  // Adjust the hit
299  if (mag(nearPoint - missPoint) < mergeTol)
300  {
301  if (debug)
302  {
303  Pout<< "a";
304  }
305 
306 // Pout<< "Moving slave point in integral adjustment "
307 // << pointi << " miss point: " << missPoint
308 // << " near point: " << nearPoint
309 // << " mergeTol: " << mergeTol
310 // << " dist: " << mag(nearPoint - missPoint) << endl;
311 
312  projectedSlavePoints[pointi] = nearPoint;
313 
314  slavePointFaceHits[pointi] =
315  objectHit(true, slavePointFaceHits[pointi].hitObject());
316 
317  nAdjustedPoints++;
318  }
319  else
320  {
321  projectedSlavePoints[pointi] = slaveLocalPoints[pointi];
322 
323  if (debug)
324  {
325  Pout<< "n";
326  }
327  }
328  }
329  }
330 
331  if (debug)
332  {
333  Pout<< " done." << endl;
334  }
335  }
336  else if (matchType_ == PARTIAL)
337  {
338  forAll(slavePointFaceHits, pointi)
339  {
340  if (slavePointFaceHits[pointi].hit())
341  {
342  // Grab the hit point
343  projectedSlavePoints[pointi] =
344  masterLocalFaces
345  [slavePointFaceHits[pointi].hitObject()].ray
346  (
347  slaveLocalPoints[pointi],
348  slavePointNormals[pointi],
349  masterLocalPoints,
350  projectionAlgo_
351  ).hitPoint();
352  }
353  else
354  {
355  // The point remains where it started from
356  projectedSlavePoints[pointi] = slaveLocalPoints[pointi];
357  }
358  }
359  }
360  else
361  {
363  << " for object " << name()
364  << abort(FatalError);
365  }
366 
367  if (debug)
368  {
369  Pout<< "Number of adjusted points in projection: "
370  << nAdjustedPoints << endl;
371 
372  // Check for zero-length edges in slave projection
373  scalar minEdgeLength = GREAT;
374  label nShortEdges = 0;
375 
376  for (const edge& e : slaveEdges)
377  {
378  const scalar len = e.mag(projectedSlavePoints);
379  minEdgeLength = min(minEdgeLength, len);
380 
381  if (len < SMALL)
382  {
383  Pout<< "Point projection problems for edge: "
384  << e << ". Length = " << len
385  << endl;
386 
387  nShortEdges++;
388  }
389  }
390 
391  if (nShortEdges > 0)
392  {
394  << nShortEdges << " short projected edges "
395  << "after adjustment for object " << name()
396  << abort(FatalError);
397  }
398  else
399  {
400  Pout<< " ... projection OK." << endl;
401  }
402  }
403 // scalarField magDiffs(mag(slaveLocalPoints - projectedSlavePoints));
404 // Pout<< "Slave point face hits: " << slavePointFaceHits << nl
405 // << "slave points: " << slaveLocalPoints << nl
406 // << "Projected slave points: " << projectedSlavePoints
407 // << "diffs: " << magDiffs << endl;
408 
409  // Merge projected points against master points
410 
411  labelList slavePointPointHits(slaveLocalPoints.size(), -1);
412  labelList masterPointPointHits(masterLocalPoints.size(), -1);
413 
414  // Go through all the slave points and compare them against all the points
415  // in the master face they hit. If the distance between the projected point
416  // and any of the master face points is smaller than both tolerances,
417  // merge the projected point by:
418  // 1) adjusting the projected point to coincide with the
419  // master point it merges with
420  // 2) remembering the hit master point index in slavePointPointHits
421  // 3) resetting the hit face to -1
422  // 4) If a master point has been hit directly, it cannot be merged
423  // into the edge. Mark it as used in the reverse list
424 
425  label nMergedPoints = 0;
426 
427  forAll(projectedSlavePoints, pointi)
428  {
429  if (slavePointFaceHits[pointi].hit())
430  {
431  // Use non-const reference so the point can be adjusted
432  point& curPoint = projectedSlavePoints[pointi];
433 
434  // Get the hit face (on master)
435  const face& hitFace =
436  masterLocalFaces[slavePointFaceHits[pointi].hitObject()];
437 
438  label mergePoint = -1;
439  scalar mergeDist = GREAT;
440 
441  // Try all point before deciding on best fit.
442  for (const label hitPointi : hitFace)
443  {
444  const scalar dist =
445  mag(masterLocalPoints[hitPointi] - curPoint);
446 
447  // Calculate the tolerance
448  const scalar mergeTol =
449  pointMergeTol_*
450  min
451  (
452  minSlavePointLength[pointi],
453  minMasterPointLength[hitPointi]
454  );
455 
456  if (dist < mergeTol && dist < mergeDist)
457  {
458  mergePoint = hitPointi;
459  mergeDist = dist;
460 
461 // Pout<< "Merging slave point "
462 // << slavePatch.meshPoints()[pointi] << " at "
463 // << slaveLocalPoints[pointi] << " with master "
464 // << masterPatch.meshPoints()[mergePoint] << " at "
465 // << masterLocalPoints[mergePoint]
466 // << ". dist: " << mergeDist
467 // << " mergeTol: " << mergeTol << endl;
468  }
469  }
470 
471  if (mergePoint > -1)
472  {
473  // Slave point is to be merged with master point.
474  // This may also include a false positive when the two points
475  // already point to the same global point, but this will need
476  // to be addressed by the caller.
477  nMergedPoints++;
478 
479  slavePointPointHits[pointi] = mergePoint;
480  masterPointPointHits[mergePoint] = pointi;
481 
482  // Adjust (snap) slave point
483  curPoint = masterLocalPoints[mergePoint];
484  }
485  }
486  }
487 
488 // Pout<< "slavePointPointHits: " << slavePointPointHits << nl
489 // << "masterPointPointHits: " << masterPointPointHits << endl;
490 
491  if (debug)
492  {
493  // Check for zero-length edges in slave projection
494  scalar minEdgeLength = GREAT;
495 
496  for (const edge& e : slaveEdges)
497  {
498  const scalar len = e.mag(projectedSlavePoints);
499  minEdgeLength = min(minEdgeLength, len);
500 
501  if (len < SMALL)
502  {
503  Pout<< "Point projection problems for edge: "
504  << e << ". Length = " << len
505  << endl;
506  }
507  }
508 
509  if (minEdgeLength < SMALL)
510  {
512  << " after point merge for object " << name()
513  << abort(FatalError);
514  }
515  else
516  {
517  Pout<< " ... point merge OK." << endl;
518  }
519  }
520 
521  // Merge projected points against master edges
522 
523  labelList slavePointEdgeHits(slaveLocalPoints.size(), -1);
524 
525  label nMovedPoints = 0;
526 
527  forAll(projectedSlavePoints, pointi)
528  {
529  // Eliminate the points merged into points
530  if (slavePointPointHits[pointi] < 0)
531  {
532  // Get current point position
533  point& curPoint = projectedSlavePoints[pointi];
534 
535  // Get the hit face
536  const labelList& hitFaceEdges =
537  masterFaceEdges[slavePointFaceHits[pointi].hitObject()];
538 
539  // Calculate the tolerance
540  const scalar mergeLength =
541  min
542  (
543  minSlavePointLength[pointi],
544  minMasterFaceLength[slavePointFaceHits[pointi].hitObject()]
545  );
546 
547  const scalar mergeTol = pointMergeTol_*mergeLength;
548 
549  scalar minDistance = GREAT;
550 
551  for (const label edgei : hitFaceEdges)
552  {
553  const edge& curEdge = masterEdges[edgei];
554 
555  pointHit edgeHit =
556  curEdge.line(masterLocalPoints).nearestDist(curPoint);
557 
558  if (edgeHit.hit())
559  {
560  const scalar dist =
561  mag(edgeHit.hitPoint() - projectedSlavePoints[pointi]);
562 
563  if (dist < mergeTol && dist < minDistance)
564  {
565  // Point is to be moved onto master edge
566  nMovedPoints++;
567 
568  slavePointEdgeHits[pointi] = edgei;
569  projectedSlavePoints[pointi] = edgeHit.hitPoint();
570 
571  minDistance = dist;
572 
573 // Pout<< "Moving slave point "
574 // << slavePatch.meshPoints()[pointi]
575 // << " (" << pointi
576 // << ") at " << slaveLocalPoints[pointi]
577 // << " onto master edge " << edgei
578 // << " or ("
579 // << masterLocalPoints[curEdge.start()]
580 // << masterLocalPoints[curEdge.end()]
581 // << ") hit: " << edgeHit.hitPoint()
582 // << ". dist: " << dist
583 // << " mergeTol: " << mergeTol << endl;
584  }
585  }
586  } // end of hit face edges
587 
588  if (slavePointEdgeHits[pointi] > -1)
589  {
590  // Projected slave point has moved. Re-attempt merge with
591  // master points of the edge
592  point& curPoint = projectedSlavePoints[pointi];
593 
594  const edge& hitMasterEdge =
595  masterEdges[slavePointEdgeHits[pointi]];
596 
597  label mergePoint = -1;
598  scalar mergeDist = GREAT;
599 
600  forAll(hitMasterEdge, hmeI)
601  {
602  const scalar hmeDist =
603  mag(masterLocalPoints[hitMasterEdge[hmeI]] - curPoint);
604 
605  // Calculate the tolerance
606  const scalar mergeTol =
607  pointMergeTol_*
608  min
609  (
610  minSlavePointLength[pointi],
611  minMasterPointLength[hitMasterEdge[hmeI]]
612  );
613 
614  if (hmeDist < mergeTol && hmeDist < mergeDist)
615  {
616  mergePoint = hitMasterEdge[hmeI];
617  mergeDist = hmeDist;
618 
619 // Pout<< "Merging slave point; SECOND TRY "
620 // << slavePatch.meshPoints()[pointi] << " local "
621 // << pointi << " at "
622 // << slaveLocalPoints[pointi] << " with master "
623 // << masterPatch.meshPoints()[mergePoint] << " at "
624 // << masterLocalPoints[mergePoint]
625 // << ". hmeDist: " << mergeDist
626 // << " mergeTol: " << mergeTol << endl;
627  }
628  }
629 
630  if (mergePoint > -1)
631  {
632  // Point is to be merged with master point
633  slavePointPointHits[pointi] = mergePoint;
634  curPoint = masterLocalPoints[mergePoint];
635  masterPointPointHits[mergePoint] = pointi;
636 
637  slavePointFaceHits[pointi] =
638  objectHit(true, slavePointFaceHits[pointi].hitObject());
639 
640 
641  // Disable edge merge
642  slavePointEdgeHits[pointi] = -1;
643  }
644  }
645  }
646  }
647 
648  if (debug)
649  {
650  Pout<< "Number of merged master points: " << nMergedPoints << nl
651  << "Number of adjusted slave points: " << nMovedPoints << endl;
652 
653  // Check for zero-length edges in slave projection
654  scalar minEdgeLength = GREAT;
655 
656  for (const edge& e : slaveEdges)
657  {
658  const scalar len = e.mag(projectedSlavePoints);
659  minEdgeLength = min(minEdgeLength, len);
660 
661  if (len < SMALL)
662  {
663  Pout<< "Point projection problems for edge: "
664  << e << ". Length = " << len
665  << endl;
666  }
667  }
668 
669  if (minEdgeLength < SMALL)
670  {
672  << " after correction for object " << name()
673  << abort(FatalError);
674  }
675  }
676 
677 // Pout<< "slavePointEdgeHits: " << slavePointEdgeHits << endl;
678 
679  // Insert the master points into closest slave edge if appropriate
680 
681  // Algorithm:
682  // The face potentially interacts with all the points of the
683  // faces covering the path between its two ends. Since we are
684  // examining an arbitrary shadow of the edge on a non-Euclidian
685  // surface, it is typically quite hard to do a geometric point
686  // search (under a shadow of a straight line). Therefore, the
687  // search will be done topologically.
688  //
689  // I) Point collection
690  // -------------------
691  // 1) Grab the master faces to which the end points of the edge
692  // have been projected.
693  // 2) Starting from the face containing the edge start, grab all
694  // its points and put them into the point lookup map. Put the
695  // face onto the face lookup map.
696  // 3) If the face of the end point is on the face lookup, complete
697  // the point collection step (this is checked during insertion.
698  // 4) Start a new round of insertion. Visit all faces in the face
699  // lookup and add their neighbours which are not already on the
700  // map. Every time the new neighbour is found, add its points to
701  // the point lookup. If the face of the end point is inserted,
702  // continue with the current roundof insertion and stop the
703  // algorithm.
704  //
705  // II) Point check
706  // ---------------
707  // Grab all the points from the point collection and check them
708  // against the current edge. If the edge-to-point distance is
709  // smaller than the smallest tolerance in the game (min of
710  // master point tolerance and left and right slave face around
711  // the edge tolerance) and the nearest point hits the edge, the
712  // master point will break the slave edge. Check the actual
713  // distance of the original master position from the edge. If
714  // the distance is smaller than a fraction of slave edge
715  // length, the hit is considered valid. Store the slave edge
716  // index for every master point.
717 
718  labelList masterPointEdgeHits(masterLocalPoints.size(), -1);
719  scalarField masterPointEdgeDist(masterLocalPoints.size(), GREAT);
720 
721  // Note. "Processing slave edges" code is repeated twice in the
722  // sliding intergace coupling in order to allow the point
723  // projection to be done separately from the actual cutting.
724  // Please change consistently with coupleSlidingInterface.C
725  //
726 
727  if (debug)
728  {
729  Pout<< "Processing slave edges " << endl;
730  }
731 
732  // Create a map of faces the edge can interact with
733  labelHashSet curFaceMap
734  (
735  nFacesPerSlaveEdge_*primitiveMesh::edgesPerFace_
736  );
737 
739 
740  forAll(slaveEdges, edgei)
741  {
742  const edge& curEdge = slaveEdges[edgei];
743 
744  {
745  // Grab the faces for start and end points
746  const label startFace =
747  slavePointFaceHits[curEdge.start()].hitObject();
748 
749  const label endFace =
750  slavePointFaceHits[curEdge.end()].hitObject();
751 
752  // Pout<< "Doing edge " << edgei
753  // << " or " << curEdge
754  // << " start: "
755  // << slavePointFaceHits[curEdge.start()].hitObject()
756  // << " end "
757  // << slavePointFaceHits[curEdge.end()].hitObject()
758  // << endl;
759 
760  // If the end face is on the list, the face collection is finished
761  bool completed = false;
762 
763  if (!completed)
764  {
765  // Forward sweep
766 
767  // Clear the maps
768  curFaceMap.clear();
769  addedFaces.clear();
770 
771  // Insert the start face into the list
772  curFaceMap.insert(startFace);
773  addedFaces.insert(startFace);
774 
775  for
776  (
777  label nSweeps = 0;
778  nSweeps < edgeFaceEscapeLimit_;
779  ++nSweeps
780  )
781  {
782  completed = addedFaces.found(endFace);
783 
784  // Add all face neighbours of face in the map
785  const labelList cf(addedFaces.toc());
786  addedFaces.clear();
787 
788  for (const label cfi : cf)
789  {
790  const labelList& curNbrs = masterFaceFaces[cfi];
791 
792  for (const label nbri : curNbrs)
793  {
794  if (curFaceMap.insert(nbri))
795  {
796  addedFaces.insert(nbri);
797  }
798  }
799  }
800 
801  if (completed) break;
802 
803  if (debug)
804  {
805  Pout<< ".";
806  }
807  }
808  }
809 
810  if (!completed)
811  {
812  // Reverse sweep
813 
814  if (debug)
815  {
816  Pout<< "x";
817  }
818 
819  // It is impossible to reach the end from the start, probably
820  // due to disconnected domain. Do search in opposite direction
821 
822  addedFaces.clear();
823  curFaceMap.insert(endFace);
824  addedFaces.insert(endFace);
825 
826  for
827  (
828  label nSweeps = 0;
829  nSweeps < edgeFaceEscapeLimit_;
830  ++nSweeps
831  )
832  {
833  completed = addedFaces.found(startFace);
834 
835  // Add all face neighbours of face in the map
836  const labelList cf(addedFaces.toc());
837  addedFaces.clear();
838 
839  for (const label cfi : cf)
840  {
841  const labelList& curNbrs = masterFaceFaces[cfi];
842 
843  for (const label nbri : curNbrs)
844  {
845  if (curFaceMap.insert(nbri))
846  {
847  addedFaces.insert(nbri);
848  }
849  }
850  }
851 
852  if (completed) break;
853 
854  if (debug)
855  {
856  Pout<< ".";
857  }
858  }
859  }
860 
861  if (completed)
862  {
863  if (debug)
864  {
865  Pout<< "+ ";
866  }
867  }
868  else
869  {
870  if (debug)
871  {
872  Pout<< "z ";
873  }
874  }
875 
876  // Collect the points
877 
878  // Create a map of points the edge can interact with
879  labelHashSet curPointMap
880  (
881  nFacesPerSlaveEdge_*primitiveMesh::pointsPerFace_
882  );
883 
884  for (const label facei : curFaceMap)
885  {
886  const face& f = masterLocalFaces[facei];
887  curPointMap.insert(f); // Insert all face points
888  }
889 
890  const labelList curMasterPoints = curPointMap.toc();
891 
892  // Check all the points against the edge.
893 
894  linePointRef edgeLine = curEdge.line(projectedSlavePoints);
895 
896  const vector edgeVec = edgeLine.vec();
897  const scalar edgeMag = edgeLine.mag();
898 
899  // Calculate actual distance involved in projection. This
900  // is used to reject master points out of reach.
901  // Calculated as a combination of travel distance in projection and
902  // edge length
903  const scalar slaveCatchDist =
904  edgeMasterCatchFraction_*edgeMag
905  + 0.5*
906  (
907  mag
908  (
909  projectedSlavePoints[curEdge.start()]
910  - slaveLocalPoints[curEdge.start()]
911  )
912  + mag
913  (
914  projectedSlavePoints[curEdge.end()]
915  - slaveLocalPoints[curEdge.end()]
916  )
917  );
918 
919  // The point merge distance needs to be measured in the
920  // plane of the slave edge. The unit vector is calculated
921  // as a cross product of the edge vector and the edge
922  // projection direction. When checking for the distance
923  // in plane, a minimum of the master-to-edge and
924  // projected-master-to-edge distance is used, to avoid
925  // problems with badly defined master planes. HJ,
926  // 17/Oct/2004
927  const vector edgeNormalInPlane =
928  normalised
929  (
930  edgeVec
931  ^ (
932  slavePointNormals[curEdge.start()]
933  + slavePointNormals[curEdge.end()]
934  )
935  );
936 
937  for (const label cmp : curMasterPoints)
938  {
939  // Skip the current point if the edge start or end has
940  // been adjusted onto in
941  if
942  (
943  slavePointPointHits[curEdge.start()] == cmp
944  || slavePointPointHits[curEdge.end()] == cmp
945  || masterPointPointHits[cmp] > -1
946  )
947  {
948 // Pout<< "Edge already snapped to point. Skipping." << endl;
949  continue;
950  }
951 
952  // Check if the point actually hits the edge within bounds
953  pointHit edgeLineHit =
954  edgeLine.nearestDist(masterLocalPoints[cmp]);
955 
956  if (edgeLineHit.hit())
957  {
958  // If the distance to the line is smaller than
959  // the tolerance the master point needs to be
960  // inserted into the edge
961 
962  // Strict checking of slave cut to avoid capturing
963  // end points.
964  const scalar cutOnSlave =
965  ((edgeLineHit.hitPoint() - edgeLine.start()) & edgeVec)
966  /sqr(edgeMag);
967 
968  const scalar distInEdgePlane =
969  min
970  (
971  edgeLineHit.distance(),
972  mag
973  (
974  (
975  masterLocalPoints[cmp]
976  - edgeLineHit.hitPoint()
977  )
978  & edgeNormalInPlane
979  )
980  );
981 // Pout<< "master point: " << cmp
982 // << " cutOnSlave " << cutOnSlave
983 // << " distInEdgePlane: " << distInEdgePlane
984 // << " tol1: " << pointMergeTol_*edgeMag
985 // << " hitDist: " << edgeLineHit.distance()
986 // << " tol2: " <<
987 // min
988 // (
989 // slaveCatchDist,
990 // masterPointEdgeDist[cmp]
991 // ) << endl;
992 
993  // Not a point hit, check for edge
994  if
995  (
996  cutOnSlave > edgeEndCutoffTol_
997  && cutOnSlave < 1.0 - edgeEndCutoffTol_ // check edge cut
998  && distInEdgePlane < edgeMergeTol_*edgeMag // merge plane
999  && edgeLineHit.distance()
1000  < min
1001  (
1002  slaveCatchDist,
1003  masterPointEdgeDist[cmp]
1004  )
1005  )
1006  {
1007  if (debug)
1008  {
1009  if (masterPointEdgeHits[cmp] == -1)
1010  {
1011  // First hit
1012  Pout<< "m";
1013  }
1014  else
1015  {
1016  // Repeat hit
1017  Pout<< "M";
1018  }
1019  }
1020 
1021  // Snap to point onto edge
1022  masterPointEdgeHits[cmp] = edgei;
1023  masterPointEdgeDist[cmp] = edgeLineHit.distance();
1024 
1025 // Pout<< "Inserting master point "
1026 // << masterPatch.meshPoints()[cmp]
1027 // << " (" << cmp
1028 // << ") at " << masterLocalPoints[cmp]
1029 // << " into slave edge " << edgei
1030 // << " " << curEdge
1031 // << " cutOnSlave: " << cutOnSlave
1032 // << " distInEdgePlane: " << distInEdgePlane
1033 // << ". dist: " << edgeLineHit.distance()
1034 // << " mergeTol: "
1035 // << edgeMergeTol_*edgeMag
1036 // << " other tol: " <<
1037 // min
1038 // (
1039 // slaveCatchDist,
1040 // masterPointEdgeDist[cmp]
1041 // )
1042 // << endl;
1043  }
1044  }
1045  }
1046  } // End if both ends missing
1047  } // End all slave edges
1048 
1049  if (debug)
1050  {
1051  Pout<< endl;
1052  }
1053 // Pout<< "masterPointEdgeHits: " << masterPointEdgeHits << endl;
1054 
1055  if (debug)
1056  {
1057  Pout<< "bool slidingInterface::projectPoints() for object "
1058  << name() << " : "
1059  << "Finished projecting points. Topology = ";
1060  }
1061 
1062 // Pout<< "slavePointPointHits: " << slavePointPointHits << nl
1063 // << "slavePointEdgeHits: " << slavePointEdgeHits << nl
1064 // << "slavePointFaceHits: " << slavePointFaceHits << nl
1065 // << "masterPointEdgeHits: " << masterPointEdgeHits << endl;
1066 
1067  // The four lists:
1068  // - slavePointPointHits
1069  // - slavePointEdgeHits
1070  // - slavePointFaceHits
1071  // - masterPointEdgeHits
1072  // define how the two patches will be merged topologically.
1073  // If the lists have not changed since the last merge, the
1074  // sliding interface changes only geometrically and simple mesh
1075  // motion will suffice. Otherwise, a topological change is
1076  // required.
1077 
1078  // Compare the result with the current state
1079  if (!attached_)
1080  {
1081  // The mesh needs to change topologically
1082  trigger_ = true;
1083 
1084  // Store the addressing arrays and projected points
1085  slavePointPointHitsPtr_.reset
1086  (
1087  new labelList(std::move(slavePointPointHits))
1088  );
1089  slavePointEdgeHitsPtr_.reset
1090  (
1091  new labelList(std::move(slavePointEdgeHits))
1092  );
1093  slavePointFaceHitsPtr_.reset
1094  (
1095  new List<objectHit>(std::move(slavePointFaceHits))
1096  );
1097  masterPointEdgeHitsPtr_.reset
1098  (
1099  new labelList(std::move(masterPointEdgeHits))
1100  );
1101 
1102  if (debug)
1103  {
1104  Pout<< "(Detached interface) changing." << endl;
1105  }
1106  }
1107  else
1108  {
1109  // Compare the lists against the stored lists. If any of them
1110  // has changed, topological change will be executed.
1111  trigger_ = false;
1112 
1113  if
1114  (
1115  !slavePointPointHitsPtr_
1116  || !slavePointEdgeHitsPtr_
1117  || !slavePointFaceHitsPtr_
1118  || !masterPointEdgeHitsPtr_
1119  )
1120  {
1121  // Must be restart. Force topological change.
1122  slavePointPointHitsPtr_.reset
1123  (
1124  new labelList(slavePointPointHits)
1125  );
1126  slavePointEdgeHitsPtr_.reset
1127  (
1128  new labelList(std::move(slavePointEdgeHits))
1129  );
1130  slavePointFaceHitsPtr_.reset
1131  (
1132  new List<objectHit>(std::move(slavePointFaceHits))
1133  );
1134  masterPointEdgeHitsPtr_.reset
1135  (
1136  new labelList(std::move(masterPointEdgeHits))
1137  );
1138 
1139  if (debug)
1140  {
1141  Pout<< "(Attached interface restart) changing." << endl;
1142  }
1143 
1144  trigger_ = true;
1145  return trigger_;
1146  }
1147 
1148  if (slavePointPointHits != *slavePointPointHitsPtr_)
1149  {
1150  if (debug)
1151  {
1152  Pout<< "(Point projection) ";
1153  }
1154 
1155  trigger_ = true;
1156  }
1157 
1158  if (slavePointEdgeHits != *slavePointEdgeHitsPtr_)
1159  {
1160  if (debug)
1161  {
1162  Pout<< "(Edge projection) ";
1163  }
1164 
1165  trigger_ = true;
1166  }
1167 
1168  // Comparison for face hits needs to exclude points that have hit
1169  // another point or edge
1170  bool faceHitsDifferent = false;
1171 
1172  const List<objectHit>& oldPointFaceHits = *slavePointFaceHitsPtr_;
1173 
1174  forAll(slavePointFaceHits, pointi)
1175  {
1176  if
1177  (
1178  slavePointPointHits[pointi] < 0
1179  && slavePointEdgeHits[pointi] < 0
1180  )
1181  {
1182  // This is a straight face hit
1183  if (slavePointFaceHits[pointi] != oldPointFaceHits[pointi])
1184  {
1185  // Two lists are different
1186  faceHitsDifferent = true;
1187  break;
1188  }
1189  }
1190  }
1191 
1192  if (faceHitsDifferent)
1193  {
1194  if (debug)
1195  {
1196  Pout<< "(Face projection) ";
1197  }
1198 
1199  trigger_ = true;
1200 
1201  }
1202 
1203  if (masterPointEdgeHits != *masterPointEdgeHitsPtr_)
1204  {
1205  if (debug)
1206  {
1207  Pout<< "(Master point projection) ";
1208  }
1209 
1210  trigger_ = true;
1211  }
1212 
1213 
1214  if (trigger_)
1215  {
1216  // Grab new data
1217  slavePointPointHitsPtr_.reset
1218  (
1219  new labelList(std::move(slavePointPointHits))
1220  );
1221  slavePointEdgeHitsPtr_.reset
1222  (
1223  new labelList(std::move(slavePointEdgeHits))
1224  );
1225  slavePointFaceHitsPtr_.reset
1226  (
1227  new List<objectHit>(std::move(slavePointFaceHits))
1228  );
1229  masterPointEdgeHitsPtr_.reset
1230  (
1231  new labelList(std::move(masterPointEdgeHits))
1232  );
1233 
1234  if (debug)
1235  {
1236  Pout<< "changing." << endl;
1237  }
1238  }
1239  else
1240  {
1241  if (debug)
1242  {
1243  Pout<< "preserved." << endl;
1244  }
1245  }
1246  }
1247 
1248  return trigger_;
1249 }
1250 
1251 
1252 void Foam::slidingInterface::clearPointProjection() const
1253 {
1254  slavePointPointHitsPtr_.reset(nullptr);
1255  slavePointEdgeHitsPtr_.reset(nullptr);
1256  slavePointFaceHitsPtr_.reset(nullptr);
1257  masterPointEdgeHitsPtr_.reset(nullptr);
1258 
1259  projectedSlavePointsPtr_.reset(nullptr);
1260 }
1261 
1262 
1263 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::primitiveMesh::pointsPerFace_
static const unsigned pointsPerFace_
Estimated number of points per face.
Definition: primitiveMesh.H:425
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::primitiveFacePatch
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
Definition: primitiveFacePatch.H:51
Foam::polyMeshModifier::topoChanger
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
Definition: polyMeshModifier.C:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::slidingInterface::PARTIAL
Definition: slidingInterface.H:86
polyTopoChanger.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::primitiveMesh::edgesPerFace_
static const unsigned edgesPerFace_
Estimated number of edges per cell.
Definition: primitiveMesh.H:416
line.H
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::FatalError
error FatalError
Foam::pointHit
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::polyMeshModifier::name
const word & name() const
Return name of this modifier.
Definition: polyMeshModifier.H:150
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
slidingInterface.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::polyTopoChanger::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyTopoChanger.H:111
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::slidingInterface::INTEGRAL
Definition: slidingInterface.H:85
Foam::DynamicID::index
label index() const
The index of the first matching items, -1 if no matches.
Definition: DynamicID.H:123