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