edgeIntersections.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) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
25 
26 #include "edgeIntersections.H"
27 #include "triSurfaceSearch.H"
28 #include "labelPairLookup.H"
29 #include "OFstream.H"
30 #include "HashSet.H"
31 #include "triSurface.H"
32 #include "pointIndexHit.H"
33 #include "treeDataTriSurface.H"
34 #include "indexedOctree.H"
35 #include "meshTools.H"
36 #include "plane.H"
37 #include "Random.H"
38 #include "unitConversion.H"
39 #include "treeBoundBox.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 defineTypeNameAndDebug(edgeIntersections, 0);
46 
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
54 {
55  const pointField& localPoints = surf.localPoints();
56  const edgeList& edges = surf.edges();
57 
58  treeBoundBox bb(localPoints);
59 
60  scalar minSize = SMALL * bb.minDim();
61 
62  forAll(edges, edgeI)
63  {
64  const edge& e = edges[edgeI];
65 
66  scalar eMag = e.mag(localPoints);
67 
68  if (eMag < minSize)
69  {
71  << "Edge " << edgeI << " vertices " << e
72  << " coords:" << localPoints[e[0]] << ' '
73  << localPoints[e[1]] << " is very small compared to bounding"
74  << " box dimensions " << bb << endl
75  << "This might lead to problems in intersection"
76  << endl;
77  }
78  }
79 }
80 
81 
82 // Update intersections for selected edges.
84 (
85  const triSurface& surf1,
86  const pointField& points1, // surf1 meshPoints (not localPoints!)
87  const triSurfaceSearch& querySurf2,
88  const scalarField& surf1PointTol, // surf1 tolerance per point
89  const labelList& edgeLabels
90 )
91 {
92  const triSurface& surf2 = querySurf2.surface();
93  const vectorField& normals2 = surf2.faceNormals();
94 
95  const labelList& meshPoints = surf1.meshPoints();
96 
97  if (debug)
98  {
99  Pout<< "Calculating intersection of " << edgeLabels.size() << " edges"
100  << " out of " << surf1.nEdges() << " with " << surf2.size()
101  << " triangles ..." << endl;
102  }
103 
104  pointField start(edgeLabels.size());
105  pointField end(edgeLabels.size());
106  vectorField edgeDirs(edgeLabels.size());
107 
108  // Go through all edges, calculate intersections
109  forAll(edgeLabels, i)
110  {
111  label edgeI = edgeLabels[i];
112 
113  if (debug)// && (i % 1000 == 0))
114  {
115  Pout<< "Intersecting edge " << edgeI << " with surface" << endl;
116  }
117 
118  const edge& e = surf1.edges()[edgeI];
119 
120  const point& pStart = points1[meshPoints[e.start()]];
121  const point& pEnd = points1[meshPoints[e.end()]];
122 
123  const vector eVec(pEnd - pStart);
124  const vector n(eVec/(mag(eVec) + VSMALL));
125 
126  // Start tracking somewhat before pStart and up to somewhat after p1.
127  // Note that tolerances here are smaller than those used to classify
128  // hit below.
129  // This will cause this hit to be marked as degenerate and resolved
130  // later on.
131  start[i] = pStart - 0.5*surf1PointTol[e[0]]*n;
132  end[i] = pEnd + 0.5*surf1PointTol[e[1]]*n;
133 
134  edgeDirs[i] = n;
135  }
136 
138  querySurf2.findLineAll
139  (
140  start,
141  end,
143  );
144 
145  label nHits = 0;
146 
147  // Classify the hits
148  forAll(edgeLabels, i)
149  {
150  const label edgeI = edgeLabels[i];
151 
152  labelList& intersectionTypes = classification_[edgeI];
153  intersectionTypes.setSize(edgeIntersections[i].size(), -1);
154 
155  this->operator[](edgeI).transfer(edgeIntersections[i]);
156 
157  forAll(intersectionTypes, hitI)
158  {
159  const pointIndexHit& pHit = this->operator[](edgeI)[hitI];
160  label& hitType = intersectionTypes[hitI];
161 
162  if (!pHit.hit())
163  {
164  continue;
165  }
166 
167  const edge& e = surf1.edges()[edgeI];
168 
169  // Classify point on surface1 edge.
170  if (mag(pHit.hitPoint() - start[i]) < surf1PointTol[e[0]])
171  {
172  // Intersection is close to edge start
173  hitType = 0;
174  }
175  else if (mag(pHit.hitPoint() - end[i]) < surf1PointTol[e[1]])
176  {
177  // Intersection is close to edge end
178  hitType = 1;
179  }
180  else if (mag(edgeDirs[i] & normals2[pHit.index()]) < alignedCos_)
181  {
182  // Edge is almost coplanar with the face
183 
184  Pout<< "Flat angle edge:" << edgeI
185  << " face:" << pHit.index()
186  << " cos:" << mag(edgeDirs[i] & normals2[pHit.index()])
187  << endl;
188 
189  hitType = 2;
190  }
191 
192  if (debug)
193  {
194  Info<< " hit " << pHit << " classify = " << hitType << endl;
195  }
196 
197  nHits++;
198  }
199  }
200 
201  if (debug)
202  {
203  Pout<< "Found " << nHits << " intersections of edges with surface ..."
204  << endl;
205  }
206 }
207 
208 
209 // If edgeI intersections are close to endpoint of edge shift endpoints
210 // slightly along edge
211 // Updates
212 // - points1 with new endpoint position
213 // - affectedEdges with all edges affected by moving the point
214 // Returns true if changed anything.
216 (
217  const triSurface& surf1,
218  const scalarField& surf1PointTol, // surf1 tolerance per point
219  const label edgeI,
220  Random& rndGen,
221  pointField& points1,
222  boolList& affectedEdges
223 ) const
224 {
225  bool hasPerturbed = false;
226 
227  // Check if edge close to endpoint. Note that we only have to check
228  // the intersection closest to the edge endpoints (i.e. first and last in
229  // edgeEdgs)
230 
231  const labelList& edgeEnds = classification_[edgeI];
232 
233  if (edgeEnds.size())
234  {
235  bool perturbStart = false;
236  bool perturbEnd = false;
237 
238  // Check first intersection.
239  if (edgeEnds.first() == 0)
240  {
241  perturbStart = true;
242  }
243 
244  if (edgeEnds.last() == 1)
245  {
246  perturbEnd = true;
247  }
248 
249 
250  if (perturbStart || perturbEnd)
251  {
252  const edge& e = surf1.edges()[edgeI];
253 
254  label v0 = surf1.meshPoints()[e[0]];
255  label v1 = surf1.meshPoints()[e[1]];
256 
257  vector eVec(points1[v1] - points1[v0]);
258  vector n = eVec/mag(eVec);
259 
260  if (perturbStart)
261  {
262  // Perturb with something (hopefully) larger than tolerance.
263  scalar t = 4.0*(rndGen.scalar01() - 0.5);
264  points1[v0] += t*surf1PointTol[e[0]]*n;
265 
266  const labelList& pEdges = surf1.pointEdges()[e[0]];
267 
268  forAll(pEdges, i)
269  {
270  affectedEdges[pEdges[i]] = true;
271  }
272  }
273  if (perturbEnd)
274  {
275  // Perturb with something larger than tolerance.
276  scalar t = 4.0*(rndGen.scalar01() - 0.5);
277  points1[v1] += t*surf1PointTol[e[1]]*n;
278 
279  const labelList& pEdges = surf1.pointEdges()[e[1]];
280 
281  forAll(pEdges, i)
282  {
283  affectedEdges[pEdges[i]] = true;
284  }
285  }
286  hasPerturbed = true;
287  }
288  }
289 
290  return hasPerturbed;
291 }
292 
293 
294 // Perturb single edge endpoint when perpendicular to face
296 (
297  const triSurface& surf1,
298  const scalarField& surf1PointTol, // surf1 tolerance per point
299  const label edgeI,
300 
301  Random& rndGen,
302  pointField& points1,
303  boolList& affectedEdges
304 ) const
305 {
306  const labelList& meshPoints = surf1.meshPoints();
307 
308  const labelList& edgeEnds = classification_[edgeI];
309 
310  bool hasPerturbed = false;
311 
312  forAll(edgeEnds, i)
313  {
314  if (edgeEnds[i] == 2)
315  {
316  const edge& e = surf1.edges()[edgeI];
317 
318  // Endpoint to modify. Choose either start or end.
319  label pointI = e[rndGen.bit()];
320  //label pointI = e[0];
321 
322  // Generate random vector slightly larger than tolerance.
323  vector rndVec = rndGen.vector01() - vector(0.5, 0.5, 0.5);
324 
325  // Make sure rndVec only perp to edge
326  vector n(points1[meshPoints[e[1]]] - points1[meshPoints[e[0]]]);
327  scalar magN = mag(n) + VSMALL;
328  n /= magN;
329 
330  rndVec -= n*(n & rndVec);
331 
332  // Normalize
333  rndVec /= mag(rndVec) + VSMALL;
334 
335  // Scale to be moved by tolerance.
336  rndVec *= 0.01*magN;
337 
338  Pout<< "rotating: shifting endpoint " << meshPoints[pointI]
339  << " of edge:" << edgeI << " verts:"
340  << points1[meshPoints[e[0]]] << ' '
341  << points1[meshPoints[e[1]]]
342  << " by " << rndVec
343  << " tol:" << surf1PointTol[pointI] << endl;
344 
345  points1[meshPoints[pointI]] += rndVec;
346 
347  // Mark edges affected by change to point
348  const labelList& pEdges = surf1.pointEdges()[pointI];
349 
350  forAll(pEdges, i)
351  {
352  affectedEdges[pEdges[i]] = true;
353  }
354 
355  hasPerturbed = true;
356 
357  // Enough done for current edge; no need to test other intersections
358  // of this edge.
359  break;
360  }
361  }
362 
363  return hasPerturbed;
364 }
365 
366 
367 // Perturb edge when close to face
369 (
370  const triSurface& surf1,
371  const triSurface& surf2,
372  const label edgeI,
373 
374  Random& rndGen,
375  pointField& points1,
376  boolList& affectedEdges
377 ) const
378 {
379  const labelList& meshPoints = surf1.meshPoints();
380 
381  const edge& e = surf1.edges()[edgeI];
382 
383  const List<pointIndexHit>& hits = operator[](edgeI);
384 
385 
386  bool hasPerturbed = false;
387 
388  // For all hits on edge
389  forAll(hits, i)
390  {
391  const pointIndexHit& pHit = hits[i];
392 
393  // Classify point on face of surface2
394  label surf2FaceI = pHit.index();
395 
396  const triSurface::FaceType& f2 = surf2.localFaces()[surf2FaceI];
397  const pointField& surf2Pts = surf2.localPoints();
398 
399  const point ctr = f2.centre(surf2Pts);
400 
401  label nearType, nearLabel;
402 
403  f2.nearestPointClassify(pHit.hitPoint(), surf2Pts, nearType, nearLabel);
404 
405  if (nearType == triPointRef::POINT || nearType == triPointRef::EDGE)
406  {
407  // Shift edge towards tri centre
408  vector offset = 0.01*rndGen.scalar01()*(ctr - pHit.hitPoint());
409 
410  // shift e[0]
411  points1[meshPoints[e[0]]] += offset;
412 
413  // Mark edges affected by change to e0
414  const labelList& pEdges0 = surf1.pointEdges()[e[0]];
415 
416  forAll(pEdges0, i)
417  {
418  affectedEdges[pEdges0[i]] = true;
419  }
420 
421  // shift e[1]
422  points1[meshPoints[e[1]]] += offset;
423 
424  // Mark edges affected by change to e1
425  const labelList& pEdges1 = surf1.pointEdges()[e[1]];
426 
427  forAll(pEdges1, i)
428  {
429  affectedEdges[pEdges1[i]] = true;
430  }
431 
432  hasPerturbed = true;
433 
434  // No need to test any other hits on this edge
435  break;
436  }
437  }
438 
439  return hasPerturbed;
440 }
441 
442 
443 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
444 
445 // Construct null
447 :
448  List<List<pointIndexHit> >(),
449  classification_()
450 {}
451 
452 
454 (
455  const triSurface& surf1,
456  const triSurfaceSearch& query2,
457  const scalarField& surf1PointTol
458 )
459 :
460  List<List<pointIndexHit> >(surf1.nEdges()),
461  classification_(surf1.nEdges())
462 {
463  checkEdges(surf1);
464 
465  // Current set of edges to test
466  labelList edgesToTest(identity(surf1.nEdges()));
467 
468  // Determine intersections for edgesToTest
469  intersectEdges
470  (
471  surf1,
472  surf1.points(), // surf1 meshPoints (not localPoints!)
473  query2,
474  surf1PointTol,
475  edgesToTest
476  );
477 }
478 
479 
480 // Construct from components
482 (
483  const List<List<pointIndexHit> >& intersections,
484  const labelListList& classification
485 )
486 :
487  List<List<pointIndexHit> >(intersections),
488  classification_(classification)
489 {}
490 
491 
492 // * * * * * * * * * * * * * * * Static Functions * * * * * * * * * * * * * //
493 
495 {
496  const pointField& localPoints = surf.localPoints();
497  const labelListList& pointEdges = surf.pointEdges();
498  const edgeList& edges = surf.edges();
499 
500  scalarField minLen(localPoints.size());
501 
502  forAll(minLen, pointI)
503  {
504  const labelList& pEdges = pointEdges[pointI];
505 
506  scalar minDist = GREAT;
507 
508  forAll(pEdges, i)
509  {
510  minDist = min(minDist, edges[pEdges[i]].mag(localPoints));
511  }
512 
513  minLen[pointI] = minDist;
514  }
515  return minLen;
516 }
517 
518 
519 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
520 
522 (
523  const label nIters,
524  const triSurface& surf1,
525  const triSurfaceSearch& query2,
526  const scalarField& surf1PointTol,
527  pointField& points1
528 )
529 {
530  const triSurface& surf2 = query2.surface();
531 
532  Random rndGen(356574);
533 
534  // Current set of edges to (re)test
535  labelList edgesToTest(surf1.nEdges());
536 
537  // Start off with all edges
538  forAll(edgesToTest, i)
539  {
540  edgesToTest[i] = i;
541  }
542 
543 
544  label iter = 0;
545 
546  for (; iter < nIters; iter++)
547  {
548  // Go through all edges to (re)test and perturb points if they are
549  // degenerate hits. Mark off edges that need to be recalculated.
550 
551  boolList affectedEdges(surf1.nEdges(), false);
552  label nShifted = 0;
553  label nRotated = 0;
554  label nOffset = 0;
555 
556  forAll(edgesToTest, i)
557  {
558  label edgeI = edgesToTest[i];
559 
560  // If edge not already marked for retesting
561  if (!affectedEdges[edgeI])
562  {
563  // 1. Check edges close to endpoint and perturb if necessary.
564 
565  bool shiftedEdgeEndPoints =
566  inlinePerturb
567  (
568  surf1,
569  surf1PointTol,
570  edgeI,
571  rndGen,
572  points1,
573  affectedEdges
574  );
575 
576  nShifted += (shiftedEdgeEndPoints ? 1 : 0);
577 
578  if (!shiftedEdgeEndPoints)
579  {
580  bool rotatedEdge =
581  rotatePerturb
582  (
583  surf1,
584  surf1PointTol,
585  edgeI,
586  rndGen,
587  points1,
588  affectedEdges
589  );
590 
591  nRotated += (rotatedEdge ? 1 : 0);
592 
593  if (!rotatedEdge)
594  {
595  // 2. we're sure now that the edge actually pierces the
596  // face. Now check the face for intersections close its
597  // points/edges
598 
599  bool offsetEdgePoints =
600  offsetPerturb
601  (
602  surf1,
603  surf2,
604  edgeI,
605  rndGen,
606  points1,
607  affectedEdges
608  );
609 
610  nOffset += (offsetEdgePoints ? 1 : 0);
611  }
612  }
613  }
614  }
615 
616  if (debug)
617  {
618  Pout<< "Edges to test : " << nl
619  << " total:" << edgesToTest.size() << nl
620  << " resolved by:" << nl
621  << " shifting : " << nShifted << nl
622  << " rotating : " << nRotated << nl
623  << " offsetting : " << nOffset << nl
624  << endl;
625  }
626 
627 
628  if (nShifted == 0 && nRotated == 0 && nOffset == 0)
629  {
630  // Nothing changed in current iteration. Current hit pattern is
631  // without any degenerates.
632  break;
633  }
634 
635  // Repack affected edges
636  labelList newEdgesToTest(surf1.nEdges());
637  label newEdgeI = 0;
638 
639  forAll(affectedEdges, edgeI)
640  {
641  if (affectedEdges[edgeI])
642  {
643  newEdgesToTest[newEdgeI++] = edgeI;
644  }
645  }
646  newEdgesToTest.setSize(newEdgeI);
647 
648  if (debug)
649  {
650  Pout<< "Edges to test:" << nl
651  << " was : " << edgesToTest.size() << nl
652  << " is : " << newEdgesToTest.size() << nl
653  << endl;
654  }
655 
656  // Transfer and test.
657  edgesToTest.transfer(newEdgesToTest);
658 
659  if (edgesToTest.empty())
660  {
661  FatalErrorInFunction << "oops" << abort(FatalError);
662  }
663 
664  // Re intersect moved edges.
665  intersectEdges
666  (
667  surf1,
668  points1, // surf1 meshPoints (not localPoints!)
669  query2,
670  surf1PointTol,
671  edgesToTest
672  );
673  }
674 
675  return iter;
676 }
677 
678 
680 (
681  const edgeIntersections& subInfo,
682  const labelList& edgeMap,
683  const labelList& faceMap,
684  const bool merge
685 )
686 {
687  forAll(subInfo, subI)
688  {
689  const List<pointIndexHit>& subHits = subInfo[subI];
690  const labelList& subClass = subInfo.classification()[subI];
691 
692  label edgeI = edgeMap[subI];
693  List<pointIndexHit>& intersections = operator[](edgeI);
694  labelList& intersectionTypes = classification_[edgeI];
695 
696  // Count unique hits. Assume edge can hit face only once
697  label sz = 0;
698  if (merge)
699  {
700  sz = intersections.size();
701  }
702 
703  label nNew = 0;
704  forAll(subHits, i)
705  {
706  const pointIndexHit& subHit = subHits[i];
707 
708  bool foundFace = false;
709  for (label interI = 0; interI < sz; interI++)
710  {
711  if (intersections[interI].index() == faceMap[subHit.index()])
712  {
713  foundFace = true;
714  break;
715  }
716  }
717  if (!foundFace)
718  {
719  nNew++;
720  }
721  }
722 
723 
724  intersections.setSize(sz+nNew);
725  intersectionTypes.setSize(sz+nNew);
726  nNew = sz;
727 
728  forAll(subHits, i)
729  {
730  const pointIndexHit& subHit = subHits[i];
731 
732  bool foundFace = false;
733  for (label interI = 0; interI < sz; interI++)
734  {
735  if (intersections[interI].index() == faceMap[subHit.index()])
736  {
737  foundFace = true;
738  break;
739  }
740  }
741 
742  if (!foundFace)
743  {
744  intersections[nNew] = pointIndexHit
745  (
746  subHit.hit(),
747  subHit.rawPoint(),
748  faceMap[subHit.index()]
749  );
750  intersectionTypes[nNew] = subClass[i];
751  nNew++;
752  }
753  }
754  }
755 }
756 
757 
758 // ************************************************************************* //
Foam::edgeIntersections::classification
const labelListList & classification() const
For every intersection the classification status.
Definition: edgeIntersections.H:178
Foam::Random
Simple random number generator.
Definition: Random.H:49
meshTools.H
labelPairLookup.H
Foam::PointIndexHit::index
label index() const
Return index.
Definition: PointIndexHit.H:111
pointIndexHit.H
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
edgeIntersections.H
Foam::edgeIntersections::alignedCos_
static scalar alignedCos_
Cosine between edge and face normal when considered parallel.
Definition: edgeIntersections.H:142
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::triFace::centre
point centre(const pointField &) const
Return centre (centroid)
Definition: triFaceI.H:163
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
indexedOctree.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
unitConversion.H
Unit conversion functions.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::triSurfaceSearch::findLineAll
void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &info) const
Calculate all intersections from start to end.
Definition: triSurfaceSearch.C:371
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::triSurfaceSearch
Helper class to search on triSurface.
Definition: triSurfaceSearch.H:55
OFstream.H
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::edgeIntersections::inlinePerturb
bool inlinePerturb(const triSurface &surf1, const scalarField &surf1PointTol, const label edgeI, Random &rndGen, pointField &points1, boolList &affectedEdges) const
Perturb endpoints of edge if they are close to the intersection.
Definition: edgeIntersections.C:216
Foam::edgeIntersections::removeDegenerates
label removeDegenerates(const label nIters, const triSurface &surf1, const triSurfaceSearch &query2, const scalarField &surf1PointTol, pointField &points1)
Resolve ties. Shuffles points so all edge - face intersections.
Definition: edgeIntersections.C:522
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::PointIndexHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointIndexHit.H:143
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
plane.H
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
treeBoundBox.H
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::PrimitivePatch::pointEdges
const labelListList & pointEdges() const
Return point-edge addressing.
Definition: PrimitivePatchTemplate.C:332
Foam::edgeIntersections::minEdgeLength
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
Definition: edgeIntersections.C:494
Foam::boundBox::minDim
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:102
HashSet.H
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::FatalError
error FatalError
Foam::triangle::POINT
@ POINT
Definition: triangle.H:97
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::cachedRandom::scalar01
scalar scalar01()
Returns the current sample.
Definition: cachedRandom.C:32
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Random.H
Foam::edgeIntersections::intersectEdges
void intersectEdges(const triSurface &surf1, const pointField &points1, const triSurfaceSearch &querySurf2, const scalarField &surf1PointTol, const labelList &edgeLabels)
Intersect selected surface edges (edgeLabels) with surface2.
Definition: edgeIntersections.C:84
Foam::edgeIntersections::edgeIntersections
edgeIntersections()
Construct null.
Definition: edgeIntersections.C:446
Foam::List::setSize
void setSize(const label)
Reset size of List.
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Foam::triangle::EDGE
@ EDGE
Definition: triangle.H:98
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::triSurfaceSearch::surface
const triSurface & surface() const
Return reference to the surface.
Definition: triSurfaceSearch.H:125
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
triSurfaceSearch.H
Foam::edgeIntersections::checkEdges
static void checkEdges(const triSurface &surf)
Check for too small edges.
Definition: edgeIntersections.C:53
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::edgeIntersections::offsetPerturb
bool offsetPerturb(const triSurface &surf1, const triSurface &surf2, const label edgeI, Random &rndGen, pointField &points1, boolList &affectedEdges) const
Perturb edge by shifting in direction trianglecentre - intersection.
Definition: edgeIntersections.C:369
Foam::edgeIntersections::merge
void merge(const edgeIntersections &, const labelList &edgeMap, const labelList &faceMap, const bool merge=true)
Merge (or override) edge intersection for a subset.
Definition: edgeIntersections.C:680
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::edgeIntersections
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
Definition: edgeIntersections.H:60
Foam::triFace::nearestPointClassify
pointHit nearestPointClassify(const point &p, const pointField &points, label &nearType, label &nearLabel) const
Return nearest point to face and classify it:
Definition: triFaceI.H:301
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
rndGen
cachedRandom rndGen(label(0), -1)
merge
bool merge(dictionary &, const dictionary &, const bool, const HashTable< wordList, word > &)
Definition: changeDictionary.C:222
Foam::edgeIntersections::rotatePerturb
bool rotatePerturb(const triSurface &surf1, const scalarField &surf1PointTol, const label edgeI, Random &rndGen, pointField &points1, boolList &affectedEdges) const
Perturb single endpoint of edge if edge is algigned with face.
Definition: edgeIntersections.C:296
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256