intersectionPatchToPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "triIntersect.H"
28 #include "unitConversion.H"
29 #include "vtkWritePolyData.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace patchToPatches
37 {
40 
42  debug::debugSwitch((intersection::typeName + "SrcFace").c_str(), -1);
44  debug::debugSwitch((intersection::typeName + "TgtFace").c_str(), -1);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * Private Static Member Functions * * * * * * * * * * //
50 
51 template<class Type>
53 Foam::patchToPatches::intersection::triPointValues
54 (
55  const triFace& t,
56  const UList<Type>& values
57 )
58 {
59  FixedList<Type, 3> result;
60  forAll(t, i)
61  {
62  result[i] = values[t[i]];
63  }
64  return result;
65 }
66 
67 
68 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
69 
71 (
72  const face& srcFace,
73  const pointField& srcPoints,
74  const vectorField& srcPointNormals
75 )
76 {
77  static DynamicList<point> ps;
78 
79  ps.clear();
80 
81  const scalar l = sqrt(mag(srcFace.area(srcPoints)));
82 
83  forAll(srcFace, srcFacePointi)
84  {
85  const label srcPointi = srcFace[srcFacePointi];
86 
87  const point& p = srcPoints[srcPointi];
88  const vector& n = srcPointNormals[srcPointi];
89 
90  ps.append(p - l/2*n);
91  ps.append(p + l/2*n);
92  }
93 
94  return treeBoundBox(ps);
95 }
96 
97 
98 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99 
100 Foam::treeBoundBox Foam::patchToPatches::intersection::srcBox
101 (
102  const face& srcFace,
103  const pointField& srcPoints,
104  const vectorField& srcPointNormals
105 ) const
106 {
107  return srcBoxStatic(srcFace, srcPoints, srcPointNormals);
108 }
109 
110 
111 bool Foam::patchToPatches::intersection::intersectFaces
112 (
113  const primitiveOldTimePatch& srcPatch,
114  const vectorField& srcPointNormals,
115  const vectorField& srcPointNormals0,
116  const primitiveOldTimePatch& tgtPatch,
117  const label srcFacei,
118  const label tgtFacei
119 )
120 {
121  // Quick rejection based on bound box
122  const treeBoundBox srcFaceBox =
123  srcBox
124  (
125  srcPatch.localFaces()[srcFacei],
126  srcPatch.localPoints(),
127  srcPointNormals
128  );
129  const treeBoundBox tgtFaceBox(tgtPatch.points(), tgtPatch[tgtFacei]);
130  if (!srcFaceBox.overlaps(tgtFaceBox)) return false;
131 
132  // Construct face triangulations on demand
133  if (srcTriPoints_[srcFacei].empty())
134  {
135  triEngine_.triangulate
136  (
138  (
139  srcPatch.localPoints(),
140  srcPatch.localFaces()[srcFacei]
141  )
142  );
143 
144  srcTriPoints_[srcFacei] =
145  triEngine_.triPoints(srcPatch.localFaces()[srcFacei]);
146  srcTriFaceEdges_[srcFacei] = triEngine_.triEdges();
147  }
148  if (tgtTriPoints_[tgtFacei].empty())
149  {
150  triEngine_.triangulate
151  (
152  UIndirectList<point>
153  (
154  tgtPatch.localPoints(),
155  tgtPatch.localFaces()[tgtFacei]
156  )
157  );
158 
159  tgtTriPoints_[tgtFacei] =
160  triEngine_.triPoints(tgtPatch.localFaces()[tgtFacei]);
161  tgtTriFaceEdges_[tgtFacei] = triEngine_.triEdges();
162  }
163 
164  // Construct and initialise workspace
165  bool srcCouples = false;
166  couple srcCouple;
167  srcFaceEdgePart_.resize(srcPatch[srcFacei].size());
168  forAll(srcFaceEdgePart_, srcFaceEdgei)
169  {
170  const edge e =
171  srcPatch.localFaces()[srcFacei].faceEdge(srcFaceEdgei);
172  const vector eC = e.centre(srcPatch.localPoints());
173  srcFaceEdgePart_[srcFaceEdgei] = part(Zero, eC);
174  }
175 
176  bool tgtCouples = false;
177  couple tgtCouple;
178  tgtFaceEdgePart_.resize(tgtPatch[tgtFacei].size());
179  forAll(tgtFaceEdgePart_, tgtFaceEdgei)
180  {
181  const edge e =
182  tgtPatch.localFaces()[tgtFacei].faceEdge(tgtFaceEdgei);
183  const vector eC = e.centre(tgtPatch.localPoints());
184  tgtFaceEdgePart_[tgtFaceEdgei] = part(Zero, eC);
185  }
186 
187  part errorPart(Zero, srcPatch.faceCentres()[srcFacei]);
188 
189  // Cache the face area magnitudes
190  const scalar srcMagA = mag(srcPatch.faceAreas()[srcFacei]);
191  const scalar tgtMagA = mag(tgtPatch.faceAreas()[tgtFacei]);
192 
193  // Determine whether or not to debug this tri intersection
194  const bool debugTriIntersect =
195  (debugSrcFacei != -1 || debugTgtFacei != -1)
196  && (debugSrcFacei == -1 || debugSrcFacei == srcFacei)
197  && (debugTgtFacei == -1 || debugTgtFacei == tgtFacei);
198 
199  // Loop the face triangles and compute the intersections
200  bool anyCouples = false;
201  forAll(srcTriPoints_[srcFacei], srcFaceTrii)
202  {
203  const triFace& srcT = srcTriPoints_[srcFacei][srcFaceTrii];
204 
205  const FixedList<point, 3> srcPs =
206  triPointValues(srcT, srcPatch.localPoints());
207  const FixedList<vector, 3> srcNs =
208  triPointValues(srcT, srcPointNormals);
209 
210  forAll(tgtTriPoints_[tgtFacei], tgtFaceTrii)
211  {
212  const triFace tgtT =
213  reverse()
214  ? tgtTriPoints_[tgtFacei][tgtFaceTrii].reverseFace()
215  : tgtTriPoints_[tgtFacei][tgtFaceTrii];
216 
217  const FixedList<point, 3> tgtPs =
218  triPointValues(tgtT, tgtPatch.localPoints());
219 
220  // Do tri-intersection
221  ictSrcPoints_.clear();
222  ictSrcPointNormals_.clear();
223  ictTgtPoints_.clear();
224  ictPointLocations_.clear();
226  (
227  srcPs,
228  srcNs,
229  {false, false, false},
230  {-1, -1, -1},
231  tgtPs,
232  {false, false, false},
233  {-1, -1, -1},
234  ictSrcPoints_,
235  ictSrcPointNormals_,
236  ictTgtPoints_,
237  ictPointLocations_,
238  debugTriIntersect,
239  debugTriIntersect
240  ? word
241  (
242  typeName
243  + "_srcFace=" + Foam::name(srcFacei)
244  + "_tgtFace=" + Foam::name(tgtFacei)
245  + "_intersection=" + Foam::name
246  (srcFaceTrii*tgtTriPoints_[tgtFacei].size() + tgtFaceTrii)
247  )
248  : word::null
249  );
250 
251  // If there is no intersection then continue
252  if (ictPointLocations_.empty())
253  {
254  continue;
255  }
256 
257  // Mark that there has been an intersection
258  anyCouples = true;
259 
260  // Compute the intersection geometry
261  const part ictSrcPart(ictSrcPoints_);
262  const part ictTgtPart(ictTgtPoints_);
263 
264  // If the intersection is below tolerance then continue
265  if
266  (
267  mag(ictSrcPart.area) < small*srcMagA
268  || mag(ictTgtPart.area) < small*tgtMagA
269  )
270  {
271  continue;
272  }
273 
274  // Mark that the source and target faces intersect
275  srcCouples = tgtCouples = true;
276 
277  // Store the intersection geometry
278  srcCouple += ictSrcPart;
279  srcCouple.nbr += ictTgtPart;
280  if (reverse())
281  {
282  tgtCouple += ictTgtPart;
283  tgtCouple.nbr += ictSrcPart;
284  }
285  else
286  {
287  tgtCouple -= ictTgtPart;
288  tgtCouple.nbr -= ictSrcPart;
289  }
290 
291  // Store the intersection polygons for debugging
292  const label debugSrcPoint0 = debugPoints_.size();
293  const label debugTgtPoint0 =
294  debugPoints_.size() + ictSrcPoints_.size();
295  if (debug > 1)
296  {
297  debugPoints_.append(ictSrcPoints_);
298  debugPoints_.append(ictTgtPoints_);
299  debugFaces_.append
300  (
301  debugSrcPoint0 + identityMap(ictSrcPoints_.size())
302  );
303  debugFaceSrcFaces_.append(srcFacei);
304  debugFaceTgtFaces_.append(tgtFacei);
305  debugFaceSides_.append(1);
306  debugFaces_.append
307  (
308  debugTgtPoint0 + identityMap(ictTgtPoints_.size())
309  );
310  debugFaceSrcFaces_.append(srcFacei);
311  debugFaceTgtFaces_.append(tgtFacei);
312  debugFaceSides_.append(-1);
313  }
314 
315  // Store edge and error areas
316  forAll(ictPointLocations_, i0)
317  {
318  const label i1 = ictPointLocations_.fcIndex(i0);
319 
320  // Get the locations on each end of this edge of the
321  // intersection polygon
322  const triIntersect::location l0 = ictPointLocations_[i0];
323  const triIntersect::location l1 = ictPointLocations_[i1];
324 
325  // Get the geometry for the projection of this edge
326  const part ictEdgePart
327  (
328  FixedList<point, 4>
329  ({
330  ictSrcPoints_[i0],
331  ictSrcPoints_[i1],
332  ictTgtPoints_[i1],
333  ictTgtPoints_[i0]
334  })
335  );
336 
337  // Store the "side" of the intersection that this edge
338  // corresponds to
339  label ictEdgeSide = -labelMax;
340 
341  // If this edge corresponds to an edge of the source
342  // triangle
343  if
344  (
345  l0.isSrcNotTgtPoint()
346  || l1.isSrcNotTgtPoint()
347  || (
348  l0.isIntersection()
349  && l1.isIntersection()
350  && l0.srcEdgei() == l1.srcEdgei()
351  )
352  )
353  {
354  const label srcEi =
355  l0.isSrcPoint() ? l0.srcPointi()
356  : l1.isSrcPoint() ? (l1.srcPointi() + 2) % 3
357  : l0.srcEdgei();
358 
359  const label srcFaceEdgei =
360  srcTriFaceEdges_[srcFacei][srcFaceTrii][srcEi];
361 
362  if (srcFaceEdgei < srcPatch[srcFacei].size())
363  {
364  srcFaceEdgePart_[srcFaceEdgei] += ictEdgePart;
365  ictEdgeSide = 1;
366  }
367  else
368  {
369  errorPart += ictEdgePart;
370  ictEdgeSide = 0;
371  }
372  }
373 
374  // If this edge corresponds to an edge of the target
375  // triangle
376  else if
377  (
378  l0.isTgtNotSrcPoint()
379  || l1.isTgtNotSrcPoint()
380  || (
381  l0.isIntersection()
382  && l1.isIntersection()
383  && l0.tgtEdgei() == l1.tgtEdgei()
384  )
385  )
386  {
387  const label tgtEi =
388  l0.isTgtPoint() ? (l0.tgtPointi() + 2) % 3
389  : l1.isTgtPoint() ? l1.tgtPointi()
390  : l0.tgtEdgei();
391 
392  const label tgtFaceEdgei =
393  tgtTriFaceEdges_[tgtFacei][tgtFaceTrii]
394  [reverse() ? 2 - tgtEi : tgtEi];
395 
396  if (tgtFaceEdgei < tgtPatch[tgtFacei].size())
397  {
398  tgtFaceEdgePart_[tgtFaceEdgei] +=
399  reverse() ? ictEdgePart : -ictEdgePart;
400  ictEdgeSide = -1;
401  }
402  else
403  {
404  errorPart += ictEdgePart;
405  ictEdgeSide = 0;
406  }
407  }
408 
409  // No other location combinations should be possible for an
410  // intersection without any shared points
411  else
412  {
414  << "The intersection topology " << ictPointLocations_
415  << " between triangle #" << srcFaceTrii
416  << " of source face #" << srcFacei
417  << " and triangle #" << tgtFaceTrii
418  << " of target face #" << tgtFacei
419  << " was not recognised. This is a bug."
420  << exit(FatalError);
421  }
422 
423  // Store the projected edge quadrilateral for debugging
424  if (debug > 1)
425  {
426  debugFaces_.append
427  (
428  labelList
429  ({
430  debugSrcPoint0 + i0,
431  debugSrcPoint0 + i1,
432  debugTgtPoint0 + i1,
433  debugTgtPoint0 + i0
434  })
435  );
436  debugFaceSrcFaces_.append(srcFacei);
437  debugFaceTgtFaces_.append(tgtFacei);
438  debugFaceSides_.append(ictEdgeSide);
439  }
440  }
441  }
442  }
443 
444  // If the source face couples the target, then store the intersection
445  if (srcCouples)
446  {
447  srcLocalTgtFaces_[srcFacei].append(tgtFacei);
448  srcCouples_[srcFacei].append(srcCouple);
449  }
450 
451  // If any intersection has occurred then store the edge and error parts
452  if (anyCouples)
453  {
454  forAll(srcFaceEdgeParts_[srcFacei], srcFaceEdgei)
455  {
456  srcFaceEdgeParts_[srcFacei][srcFaceEdgei] +=
457  srcFaceEdgePart_[srcFaceEdgei];
458  }
459  srcErrorParts_[srcFacei] +=
460  reverse() ? sum(tgtFaceEdgePart_) : -sum(tgtFaceEdgePart_);
461  srcErrorParts_[srcFacei] += errorPart;
462  }
463 
464  // If the target face couples the source, then store in the intersection
465  if (tgtCouples)
466  {
467  tgtLocalSrcFaces_[tgtFacei].append(srcFacei);
468  tgtCouples_[tgtFacei].append(tgtCouple);
469  }
470 
471  return anyCouples;
472 }
473 
474 
475 void Foam::patchToPatches::intersection::initialise
476 (
477  const primitiveOldTimePatch& srcPatch,
478  const vectorField& srcPointNormals,
479  const vectorField& srcPointNormals0,
480  const primitiveOldTimePatch& tgtPatch
481 )
482 {
484  (
485  srcPatch,
486  srcPointNormals,
487  srcPointNormals0,
488  tgtPatch
489  );
490 
491  srcCouples_.resize(srcPatch.size());
492  forAll(srcLocalTgtFaces_, i)
493  {
494  srcCouples_[i].clear();
495  }
496 
497  srcEdgeParts_.resize(srcPatch.nEdges());
498  forAll(srcEdgeParts_, srcEdgei)
499  {
500  const edge& e = srcPatch.edges()[srcEdgei];
501  const point c = e.centre(srcPatch.localPoints());
502  srcEdgeParts_[srcEdgei] = part(Zero, c);
503  }
504 
505  srcErrorParts_.resize(srcPatch.size());
506  forAll(srcErrorParts_, srcFacei)
507  {
508  srcErrorParts_[srcFacei] =
509  part(Zero, srcPatch.faceCentres()[srcFacei]);
510  }
511 
512  tgtCouples_.resize(tgtPatch.size());
513  forAll(tgtLocalSrcFaces_, i)
514  {
515  tgtCouples_[i].clear();
516  }
517 
518  srcTriPoints_ = List<triFaceList>(srcPatch.size());
519  srcTriFaceEdges_ = List<List<FixedList<label, 3>>>(srcPatch.size());
520  tgtTriPoints_ = List<triFaceList>(tgtPatch.size());
521  tgtTriFaceEdges_ = List<List<FixedList<label, 3>>>(tgtPatch.size());
522 
523  srcFaceEdgeParts_.resize(srcPatch.size());
524  forAll(srcFaceEdgeParts_, srcFacei)
525  {
526  srcFaceEdgeParts_[srcFacei].resize(srcPatch[srcFacei].size());
527  forAll(srcFaceEdgeParts_[srcFacei], srcFaceEdgei)
528  {
529  const label srcEdgei =
530  srcPatch.faceEdges()[srcFacei][srcFaceEdgei];
531  srcFaceEdgeParts_[srcFacei][srcFaceEdgei] = srcEdgeParts_[srcEdgei];
532  }
533  }
534 
535  if (debug > 1)
536  {
537  debugPoints_.clear();
538  debugFaces_.clear();
539  debugFaceSrcFaces_.clear();
540  debugFaceTgtFaces_.clear();
541  debugFaceSides_.clear();
542  }
543 }
544 
545 
546 Foam::labelList Foam::patchToPatches::intersection::finaliseLocal
547 (
548  const primitiveOldTimePatch& srcPatch,
549  const vectorField& srcPointNormals,
550  const vectorField& srcPointNormals0,
551  const primitiveOldTimePatch& tgtPatch
552 )
553 {
554  const labelList newToOldLocalTgtFace =
556  (
557  srcPatch,
558  srcPointNormals,
559  srcPointNormals0,
560  tgtPatch
561  );
562 
563  tgtCouples_ = List<DynamicList<couple>>(tgtCouples_, newToOldLocalTgtFace);
564 
565  return newToOldLocalTgtFace;
566 }
567 
568 
569 void Foam::patchToPatches::intersection::rDistributeTgt
570 (
571  const primitiveOldTimePatch& tgtPatch
572 )
573 {
575 
577  (
578  tgtPatch.size(),
579  tgtMapPtr_(),
580  tgtCouples_
581  );
582 }
583 
584 
585 Foam::label Foam::patchToPatches::intersection::finalise
586 (
587  const primitiveOldTimePatch& srcPatch,
588  const vectorField& srcPointNormals,
589  const vectorField& srcPointNormals0,
590  const primitiveOldTimePatch& tgtPatch,
591  const transformer& tgtToSrc
592 )
593 {
594  const label nCouples =
596  (
597  srcPatch,
598  srcPointNormals,
599  srcPointNormals0,
600  tgtPatch,
601  tgtToSrc
602  );
603 
604  // Convert face-edge-parts to edge-parts
605  labelList srcEdgeNParts(srcEdgeParts_.size(), 0);
606  forAll(srcEdgeParts_, srcEdgei)
607  {
608  const edge& e = srcPatch.edges()[srcEdgei];
609 
610  srcEdgeParts_[srcEdgei] = part();
611 
612  forAll(srcPatch.edgeFaces()[srcEdgei], i)
613  {
614  const label srcFacei = srcPatch.edgeFaces()[srcEdgei][i];
615  const label srcFaceEdgei =
616  findIndex(srcPatch.faceEdges()[srcFacei], srcEdgei);
617 
618  const edge fe =
619  srcPatch.localFaces()[srcFacei].faceEdge(srcFaceEdgei);
620 
621  if (edge::compare(e, fe) > 0)
622  {
623  srcEdgeParts_[srcEdgei] +=
624  srcFaceEdgeParts_[srcFacei][srcFaceEdgei];
625  }
626  else
627  {
628  srcEdgeParts_[srcEdgei] -=
629  srcFaceEdgeParts_[srcFacei][srcFaceEdgei];
630  }
631 
632  srcEdgeNParts[srcEdgei] ++;
633  }
634  }
635  forAll(srcEdgeParts_, srcEdgei)
636  {
637  srcEdgeParts_[srcEdgei].area /= srcEdgeNParts[srcEdgei];
638  }
639 
640  // Add the difference between the face-edge-part and the edge-part into the
641  // face-error-parts
642  forAll(srcEdgeParts_, srcEdgei)
643  {
644  const edge& e = srcPatch.edges()[srcEdgei];
645 
646  forAll(srcPatch.edgeFaces()[srcEdgei], i)
647  {
648  const label srcFacei = srcPatch.edgeFaces()[srcEdgei][i];
649  const label srcFaceEdgei =
650  findIndex(srcPatch.faceEdges()[srcFacei], srcEdgei);
651 
652  const edge fe =
653  srcPatch.localFaces()[srcFacei].faceEdge(srcFaceEdgei);
654 
655  if (edge::compare(e, fe) > 0)
656  {
657  srcErrorParts_[srcFacei] -= srcEdgeParts_[srcEdgei];
658  }
659  else
660  {
661  srcErrorParts_[srcFacei] += srcEdgeParts_[srcEdgei];
662  }
663 
664  srcErrorParts_[srcFacei] +=
665  srcFaceEdgeParts_[srcFacei][srcFaceEdgei];
666  }
667  }
668 
669  // Transform the target couples back to the target side
670  if (!isNull(tgtToSrc))
671  {
672  forAll(tgtCouples_, tgtFacei)
673  {
674  forAll(tgtCouples_[tgtFacei], i)
675  {
676  couple& c = tgtCouples_[tgtFacei][i];
677 
678  c.area = tgtToSrc.invTransform(c.area);
679  c.centre = tgtToSrc.invTransformPosition(c.centre);
680  c.nbr.area = tgtToSrc.invTransform(c.nbr.area);
681  c.nbr.centre = tgtToSrc.invTransformPosition(c.nbr.centre);
682  }
683  }
684  }
685 
686  // Calculate coverage and total areas on both sides
687  auto coverage = []
688  (
689  const primitivePatch& patch,
690  const List<DynamicList<couple>>& couples,
691  scalar& area,
692  scalar& coupleArea,
693  List<scalar>& coverage
694  )
695  {
696  area = 0;
697  coupleArea = 0;
698  coverage.resize(patch.size());
699 
700  forAll(patch, facei)
701  {
702  const scalar magA = mag(patch.faceAreas()[facei]);
703 
704  vector aCouple = Zero;
705  forAll(couples[facei], i)
706  {
707  aCouple += couples[facei][i].area;
708  }
709  const scalar magACouple = mag(aCouple);
710 
711  area += magA;
712  coupleArea += magACouple;
713  coverage[facei] = magACouple/magA;
714  }
715 
716  reduce(area, sumOp<scalar>());
717  reduce(coupleArea, sumOp<scalar>());
718  };
719  scalar srcArea = 0, srcCoupleArea = 0;
720  scalar tgtArea = 0, tgtCoupleArea = 0;
721  coverage(srcPatch, srcCouples_, srcArea, srcCoupleArea, srcCoverage_);
722  coverage(tgtPatch, tgtCouples_, tgtArea, tgtCoupleArea, tgtCoverage_);
723 
724  // Clear the triangulation workspace
725  srcTriPoints_.clear();
726  srcTriFaceEdges_.clear();
727  tgtTriPoints_.clear();
728  tgtTriFaceEdges_.clear();
729 
730  // Clear face-edge-parts
731  srcFaceEdgePart_.clear();
732  tgtFaceEdgePart_.clear();
733  srcFaceEdgeParts_.clear();
734 
735  // Checking and reporting
736  if (nCouples != 0)
737  {
738  scalarField srcOpenness(srcPatch.size());
739  scalarField srcError(srcPatch.size());
740  scalarField srcDepth(srcPatch.size());
741  scalarField srcAngleDeg(srcPatch.size());
742  forAll(srcPatch, srcFacei)
743  {
744  const vector& a = srcPatch.faceAreas()[srcFacei];
745  const scalar magA = mag(a);
746  const point& c = srcPatch.faceCentres()[srcFacei];
747 
748  couple Cpl(part(Zero, c), part(Zero, c));
749  forAll(srcCouples_[srcFacei], srcTgtFacei)
750  {
751  const couple& cpl = srcCouples_[srcFacei][srcTgtFacei];
752 
753  Cpl += cpl;
754  Cpl.nbr += cpl.nbr;
755  }
756 
757  vector projectionA = Zero;
758  scalar projectionV = 0;
759  forAll(srcCouples_[srcFacei], srcTgtFacei)
760  {
761  const couple& cpl = srcCouples_[srcFacei][srcTgtFacei];
762 
763  projectionA += cpl.nbr.area;
764  projectionV +=
765  - (cpl.area/3 & (cpl.centre - Cpl.centre))
766  + (cpl.nbr.area/3 & (cpl.nbr.centre - Cpl.centre));
767  }
768  forAll(srcPatch.faceEdges()[srcFacei], srcFaceEdgei)
769  {
770  const label srcEdgei =
771  srcPatch.faceEdges()[srcFacei][srcFaceEdgei];
772 
773  const edge& e = srcPatch.edges()[srcEdgei];
774  const edge fe =
775  srcPatch.localFaces()[srcFacei].faceEdge(srcFaceEdgei);
776 
777  const scalar sign = edge::compare(e, fe);
778 
779  projectionA += sign*srcEdgeParts_[srcEdgei].area;
780  projectionV +=
781  sign*srcEdgeParts_[srcEdgei].area/3
782  & (srcEdgeParts_[srcEdgei].centre - Cpl.centre);
783  }
784  projectionA += srcErrorParts_[srcFacei].area;
785  projectionV +=
786  srcErrorParts_[srcFacei].area/3
787  & (srcErrorParts_[srcFacei].centre - Cpl.centre);
788 
789  const vector aHat = normalised(a);
790  const vector aOppHat = normalised(a - Cpl.area + Cpl.nbr.area);
791  srcAngleDeg[srcFacei] =
792  radToDeg(acos(min(max(aHat & aOppHat, -1), +1)));
793  srcOpenness[srcFacei] = mag(projectionA - Cpl.area)/magA;
794  srcError[srcFacei] = mag(srcErrorParts_[srcFacei].area)/magA;
795  srcDepth[srcFacei] = mag(projectionV)/pow3(sqrt(magA));
796  }
797 
798  reduce(tgtArea, sumOp<scalar>());
799  reduce(tgtCoupleArea, sumOp<scalar>());
800 
801  Info<< indent << "Source min/average/max coverage = "
802  << gMin(srcCoverage_) << '/' << srcCoupleArea/srcArea << '/'
803  << gMax(srcCoverage_) << endl
804  << indent << "Target min/average/max coverage = "
805  << gMin(tgtCoverage_) << '/' << tgtCoupleArea/tgtArea << '/'
806  << gMax(tgtCoverage_) << endl
807  << indent << "Source average openness/error/depth/angle = "
808  << gAverage(srcOpenness) << '/' << gAverage(srcError) << '/'
809  << gAverage(srcDepth) << '/' << gAverage(srcAngleDeg) << endl
810  << indent << "Source max openness/error/depth/angle = "
811  << gMax(srcOpenness) << '/' << gMax(srcError) << '/'
812  << gMax(srcDepth) << '/' << gMax(srcAngleDeg) << endl;
813 
814  if (debug)
815  {
816  word name = patchToPatch::typeName + '_' + typeName;
817 
818  if (Pstream::parRun())
819  {
820  name += "_proc" + Foam::name(Pstream::myProcNo());
821  }
822 
823  Info<< indent << "Writing source patch to "
824  << name + "_srcPatch.vtk" << endl;
826  (
827  name + "_srcPatch" + ".vtk",
828  name + "_srcPatch",
829  false,
830  srcPatch.localPoints(),
831  labelList(),
832  labelListList(),
833  srcPatch.localFaces(),
834  "coverage", false, scalarField(srcCoverage_),
835  "openness", false, srcOpenness,
836  "error", false, srcError,
837  "depth", false, srcDepth,
838  "angle", false, srcAngleDeg,
839  "normals", true, srcPointNormals
840  );
841 
842  Info<< indent << "Writing target patch to "
843  << name + "_tgtPatch.vtk" << endl;
845  (
846  name + "_tgtPatch" + ".vtk",
847  name + "_tgtPatch",
848  false,
849  tgtPatch.localPoints(),
850  labelList(),
851  labelListList(),
852  tgtPatch.localFaces(),
853  "coverage", false, scalarField(tgtCoverage_)
854  );
855 
856  if (debug > 1)
857  {
858  Info<< indent << "Writing intersected faces to "
859  << name + ".vtk" << endl;
861  (
862  name + ".vtk",
863  name,
864  false,
865  debugPoints_,
866  labelList(),
867  labelListList(),
868  debugFaces_,
869  "srcFace", false, Field<label>(debugFaceSrcFaces_),
870  "tgtFace", false, Field<label>(debugFaceTgtFaces_),
871  "side", false, Field<label>(debugFaceSides_)
872  );
873 
874  debugPoints_.clear();
875  debugFaces_.clear();
876  debugFaceSrcFaces_.clear();
877  debugFaceTgtFaces_.clear();
878  debugFaceSides_.clear();
879  }
880  }
881  }
882 
883  return nCouples;
884 }
885 
886 
888 Foam::patchToPatches::intersection::srcWeights() const
889 {
890  List<DynamicList<scalar>>* resultPtr
891  (
892  new List<DynamicList<scalar>>(srcCouples_.size())
893  );
894  List<DynamicList<scalar>>& result = *resultPtr;
895 
896  forAll(srcCouples_, srcFacei)
897  {
898  result[srcFacei].resize(srcCouples_[srcFacei].size());
899  scalar aSum = 0;
900 
901  forAll(srcCouples_[srcFacei], i)
902  {
903  const scalar a = mag(srcCouples_[srcFacei][i].area);
904  result[srcFacei][i] = a;
905  aSum += a;
906  }
907 
908  forAll(srcCouples_[srcFacei], i)
909  {
910  result[srcFacei][i] *=
911  min(max(srcCoverage_[srcFacei], small), scalar(1))/aSum;
912  }
913  }
914 
915  return tmpNrc<List<DynamicList<scalar>>>(resultPtr);
916 }
917 
918 
920 Foam::patchToPatches::intersection::tgtWeights() const
921 {
922  List<DynamicList<scalar>>* resultPtr
923  (
924  new List<DynamicList<scalar>>(tgtCouples_.size())
925  );
926  List<DynamicList<scalar>>& result = *resultPtr;
927 
928  forAll(tgtCouples_, tgtFacei)
929  {
930  result[tgtFacei].resize(tgtCouples_[tgtFacei].size());
931  scalar aSum = 0;
932 
933  forAll(tgtCouples_[tgtFacei], i)
934  {
935  const scalar a = mag(tgtCouples_[tgtFacei][i].area);
936  result[tgtFacei][i] = a;
937  aSum += a;
938  }
939 
940  forAll(tgtCouples_[tgtFacei], i)
941  {
942  result[tgtFacei][i] *=
943  min(max(tgtCoverage_[tgtFacei], small), scalar(1))/aSum;
944  }
945  }
946 
947  return tmpNrc<List<DynamicList<scalar>>>(resultPtr);
948 }
949 
950 
951 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
952 
954 :
956 
957  srcCouples_(),
958  srcCoverage_(),
959  srcEdgeParts_(),
960  srcErrorParts_(),
961  tgtCouples_(),
962  tgtCoverage_(),
963 
964  triEngine_(),
965 
966  srcTriPoints_(),
967  srcTriFaceEdges_(),
968  tgtTriPoints_(),
969  tgtTriFaceEdges_(),
970 
971  ictSrcPoints_(),
972  ictSrcPointNormals_(),
973  ictTgtPoints_(),
974  ictPointLocations_(),
975 
976  srcFaceEdgePart_(),
977  tgtFaceEdgePart_(),
978 
979  srcFaceEdgeParts_(),
980 
981  debugPoints_(),
982  debugFaces_(),
983  debugFaceSrcFaces_(),
984  debugFaceTgtFaces_(),
985  debugFaceSides_()
986 {}
987 
988 
989 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
990 
992 {}
993 
994 
995 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
const Field< PointType > & faceAreas() const
Return face areas for patch.
const Field< PointType > & points() const
Return reference to global points.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
const Field< PointType > & faceCentres() const
Return face centres for patch.
A List with indirect addressing.
Definition: UIndirectList.H:60
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static int compare(const edge &, const edge &)
Compare edges.
Definition: edgeI.H:36
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
static vector area(const PointField &ps)
Return vector area given face points.
Class to generate coupling geometry between two primitive patches.
Definition: patchToPatch.H:56
virtual void rDistributeTgt(const primitiveOldTimePatch &tgtPatch)
Send data that resulted from an intersection between the source.
Definition: patchToPatch.C:724
virtual labelList finaliseLocal(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
Finalise the intersection locally. Trims the local target patch.
Definition: patchToPatch.C:662
virtual void initialise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch)
Initialisation.
Definition: patchToPatch.C:640
virtual label finalise(const primitiveOldTimePatch &srcPatch, const vectorField &srcPointNormals, const vectorField &srcPointNormals0, const primitiveOldTimePatch &tgtPatch, const transformer &tgtToSrc)
Finalising.
Definition: patchToPatch.C:737
Class to generate patchToPatch coupling geometry. A full geometric intersection is done between a fac...
static treeBoundBox srcBoxStatic(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals)
Get the bound box for a source face.
static int debugSrcFacei
Extra debugging for intersections between specific faces. Named.
intersection(const bool reverse)
Construct from components.
A class for managing temporary objects without reference counting.
Definition: tmpNrc.H:53
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
bool overlaps(const boundBox &) const
Overlaps other bounding box?
Definition: boundBoxI.H:126
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const dimensionedScalar c
Speed of light in a vacuum.
int debugSwitch(const char *name, const int defaultValue=0)
Lookup debug switch or add default value.
Definition: debug.C:211
static void rDistributeListList(const label size, const distributionMap &map, List< List< Type >> &data)
Reverse distribute a list-list given the map.
addToRunTimeSelectionTable(patchToPatch, intersection, bool)
defineTypeNameAndDebug(intersection, 0)
void intersectTris(const FixedList< point, 3 > &srcPs, const FixedList< vector, 3 > &srcNs, const FixedList< bool, 3 > &srcOwns, const FixedList< label, 3 > &srcTgtPis, const FixedList< point, 3 > &tgtPs, const FixedList< bool, 3 > &tgtOwns, const FixedList< label, 3 > &tgtSrcPis, DynamicList< point > &srcPoints, DynamicList< vector > &srcPointNormals, DynamicList< point > &tgtPoints, DynamicList< location > &pointLocations, const bool debug, const word &writePrefix=word::null)
Construct the intersection of a source triangle's projected volume and a.
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
scalar radToDeg(const scalar rad)
Convert radians to degrees.
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedScalar sign(const dimensionedScalar &ds)
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 & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
PrimitiveOldTimePatch< SubList< face >, const pointField & > primitiveOldTimePatch
Addressing for a faceList slice.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:510
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Type gAverage(const FieldField< Field, Type > &f)
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
static const label labelMax
Definition: label.H:62
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
Type gMax(const FieldField< Field, Type > &f)
dimensionedScalar acos(const dimensionedScalar &ds)
face triFace(3)
volScalarField & p
Functions with which to perform an intersection of a pair of triangles; the source and target....