Loading...
Searching...
No Matches
faceCoupleInfo.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "faceCoupleInfo.H"
30#include "polyMesh.H"
31#include "matchPoints.H"
33#include "meshTools.H"
34#include "treeDataFace.H"
36#include "OFstream.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43 const scalar faceCoupleInfo::angleTol_ = 1e-3;
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::faceCoupleInfo::writeOBJ
50(
51 const fileName& fName,
52 const edgeList& edges,
53 const pointField& points,
54 const bool compact
55)
56{
57 OFstream str(fName);
58
59 labelList pointMap;
60
61 if (compact)
62 {
63 pointMap.resize(points.size(), -1);
64
65 label newPointi = 0;
66
67 forAll(edges, edgeI)
68 {
69 const edge& e = edges[edgeI];
70
71 forAll(e, eI)
72 {
73 const label pointi = e[eI];
74
75 if (pointMap[pointi] == -1)
76 {
77 pointMap[pointi] = newPointi++;
78
79 meshTools::writeOBJ(str, points[pointi]);
80 }
81 }
82 }
83 }
84 else
85 {
86 pointMap = identity(points.size());
87
88 forAll(points, pointi)
89 {
90 meshTools::writeOBJ(str, points[pointi]);
91 }
92 }
93
94 forAll(edges, edgeI)
95 {
96 const edge& e = edges[edgeI];
97
98 str<< "l " << pointMap[e[0]]+1 << ' ' << pointMap[e[1]]+1 << nl;
99 }
100}
101
102
103void Foam::faceCoupleInfo::writeOBJ
104(
105 const fileName& fName,
106 const pointField& points0,
107 const pointField& points1
108)
109{
110 Pout<< "Writing connections as edges to " << fName << endl;
111
112 OFstream str(fName);
113
114 label vertI = 0;
115
116 forAll(points0, i)
117 {
119 vertI++;
120 meshTools::writeOBJ(str, points1[i]);
121 vertI++;
122 str << "l " << vertI-1 << ' ' << vertI << nl;
123 }
124}
125
126
127void Foam::faceCoupleInfo::writePointsFaces() const
128{
129 const indirectPrimitivePatch& m = masterPatch();
130 const indirectPrimitivePatch& s = slavePatch();
131 const primitiveFacePatch& c = cutFaces();
132
133 // Patches
134 {
135 OFstream str("masterPatch.obj");
136 Pout<< "Writing masterPatch to " << str.name() << endl;
137 meshTools::writeOBJ(str, m.localFaces(), m.localPoints());
138 }
139 {
140 OFstream str("slavePatch.obj");
141 Pout<< "Writing slavePatch to " << str.name() << endl;
142 meshTools::writeOBJ(str, s.localFaces(), s.localPoints());
143 }
144 {
145 OFstream str("cutFaces.obj");
146 Pout<< "Writing cutFaces to " << str.name() << endl;
147 meshTools::writeOBJ(str, c.localFaces(), c.localPoints());
148 }
149
150 // Point connectivity
151 {
152 Pout<< "Writing cutToMasterPoints to cutToMasterPoints.obj" << endl;
153
155 (
156 "cutToMasterPoints.obj",
157 m.localPoints(),
158 pointField(c.localPoints(), masterToCutPoints_));
159 }
160 {
161 Pout<< "Writing cutToSlavePoints to cutToSlavePoints.obj" << endl;
162
164 (
165 "cutToSlavePoints.obj",
166 s.localPoints(),
167 pointField(c.localPoints(), slaveToCutPoints_)
168 );
169 }
170
171 // Face connectivity
172 {
173 Pout<< "Writing cutToMasterFaces to cutToMasterFaces.obj" << endl;
174
175 pointField equivMasterFaces(c.size());
176
177 forAll(cutToMasterFaces(), cutFacei)
178 {
179 label masterFacei = cutToMasterFaces()[cutFacei];
180
181 if (masterFacei != -1)
182 {
183 equivMasterFaces[cutFacei] = m[masterFacei].centre(m.points());
184 }
185 else
186 {
188 << "No master face for cut face " << cutFacei
189 << " at position " << c[cutFacei].centre(c.points())
190 << endl;
191
192 equivMasterFaces[cutFacei] = Zero;
193 }
194 }
195
197 (
198 "cutToMasterFaces.obj",
199 calcFaceCentres<List>(c, cutPoints(), 0, c.size()),
200 equivMasterFaces
201 );
202 }
203
204 {
205 Pout<< "Writing cutToSlaveFaces to cutToSlaveFaces.obj" << endl;
206
207 pointField equivSlaveFaces(c.size());
208
209 forAll(cutToSlaveFaces(), cutFacei)
210 {
211 label slaveFacei = cutToSlaveFaces()[cutFacei];
212
213 equivSlaveFaces[cutFacei] = s[slaveFacei].centre(s.points());
214 }
215
217 (
218 "cutToSlaveFaces.obj",
219 calcFaceCentres<List>(c, cutPoints(), 0, c.size()),
220 equivSlaveFaces
221 );
222 }
223
224 Pout<< endl;
225}
226
227
228void Foam::faceCoupleInfo::writeEdges
229(
230 const labelList& cutToMasterEdges,
231 const labelList& cutToSlaveEdges
232) const
233{
234 const indirectPrimitivePatch& m = masterPatch();
235 const indirectPrimitivePatch& s = slavePatch();
236 const primitiveFacePatch& c = cutFaces();
237
238 // Edge connectivity
239 {
240 OFstream str("cutToMasterEdges.obj");
241 Pout<< "Writing cutToMasterEdges to " << str.name() << endl;
242
243 label vertI = 0;
244
245 forAll(cutToMasterEdges, cutEdgeI)
246 {
247 if (cutToMasterEdges[cutEdgeI] != -1)
248 {
249 const edge& masterEdge = m.edges()[cutToMasterEdges[cutEdgeI]];
250 const edge& cutEdge = c.edges()[cutEdgeI];
251
252 meshTools::writeOBJ(str, m.localPoints()[masterEdge[0]]);
253 vertI++;
254 meshTools::writeOBJ(str, m.localPoints()[masterEdge[1]]);
255 vertI++;
256 meshTools::writeOBJ(str, c.localPoints()[cutEdge[0]]);
257 vertI++;
258 meshTools::writeOBJ(str, c.localPoints()[cutEdge[1]]);
259 vertI++;
260 str << "l " << vertI-3 << ' ' << vertI-2 << nl;
261 str << "l " << vertI-3 << ' ' << vertI-1 << nl;
262 str << "l " << vertI-3 << ' ' << vertI << nl;
263 str << "l " << vertI-2 << ' ' << vertI-1 << nl;
264 str << "l " << vertI-2 << ' ' << vertI << nl;
265 str << "l " << vertI-1 << ' ' << vertI << nl;
266 }
267 }
268 }
269 {
270 OFstream str("cutToSlaveEdges.obj");
271 Pout<< "Writing cutToSlaveEdges to " << str.name() << endl;
272
273 label vertI = 0;
274
275 labelList slaveToCut(invert(s.nEdges(), cutToSlaveEdges));
276
277 forAll(slaveToCut, edgeI)
278 {
279 if (slaveToCut[edgeI] != -1)
280 {
281 const edge& slaveEdge = s.edges()[edgeI];
282 const edge& cutEdge = c.edges()[slaveToCut[edgeI]];
283
284 meshTools::writeOBJ(str, s.localPoints()[slaveEdge[0]]);
285 vertI++;
286 meshTools::writeOBJ(str, s.localPoints()[slaveEdge[1]]);
287 vertI++;
288 meshTools::writeOBJ(str, c.localPoints()[cutEdge[0]]);
289 vertI++;
290 meshTools::writeOBJ(str, c.localPoints()[cutEdge[1]]);
291 vertI++;
292 str << "l " << vertI-3 << ' ' << vertI-2 << nl;
293 str << "l " << vertI-3 << ' ' << vertI-1 << nl;
294 str << "l " << vertI-3 << ' ' << vertI << nl;
295 str << "l " << vertI-2 << ' ' << vertI-1 << nl;
296 str << "l " << vertI-2 << ' ' << vertI << nl;
297 str << "l " << vertI-1 << ' ' << vertI << nl;
298 }
299 }
300 }
301
302 Pout<< endl;
303}
304
305
306Foam::labelList Foam::faceCoupleInfo::findMappedEdges
307(
308 const edgeList& edges,
309 const labelList& pointMap,
310 const indirectPrimitivePatch& patch
311)
312{
313 labelList toPatchEdges(edges.size());
314
315 forAll(toPatchEdges, edgeI)
316 {
317 const edge& e = edges[edgeI];
318
319 label v0 = pointMap[e[0]];
320 label v1 = pointMap[e[1]];
321
322 toPatchEdges[edgeI] =
324 (
325 patch.edges(),
326 patch.pointEdges()[v0],
327 v0,
328 v1
329 );
330 }
331 return toPatchEdges;
332}
333
334
335bool Foam::faceCoupleInfo::regionEdge
336(
337 const polyMesh& slaveMesh,
338 const label slaveEdgeI
339) const
340{
341 const labelList& eFaces = slavePatch().edgeFaces()[slaveEdgeI];
342
343 if (eFaces.size() == 1)
344 {
345 return true;
346 }
347 else
348 {
349 // Count how many different patches connected to this edge.
350
351 label patch0 = -1;
352
353 forAll(eFaces, i)
354 {
355 label facei = eFaces[i];
356
357 label meshFacei = slavePatch().addressing()[facei];
358
359 label patchi = slaveMesh.boundaryMesh().whichPatch(meshFacei);
360
361 if (patch0 == -1)
362 {
363 patch0 = patchi;
364 }
365 else if (patchi != patch0)
366 {
367 // Found two different patches connected to this edge.
368 return true;
369 }
370 }
371 }
372
373 return false;
374}
375
376
377Foam::label Foam::faceCoupleInfo::mostAlignedCutEdge
378(
379 const bool report,
380 const polyMesh& slaveMesh,
381 const bool patchDivision,
382 const labelList& cutToMasterEdges,
383 const labelList& cutToSlaveEdges,
384 const label pointi,
385 const label edgeStart,
386 const label edgeEnd
387) const
388{
389 // Find edge using pointi that is most aligned with vector between master
390 // points. Patchdivision tells us whether or not to use patch information to
391 // match edges.
392
393 const pointField& localPoints = cutFaces().localPoints();
394
395 const labelList& pEdges = cutFaces().pointEdges()[pointi];
396
397 if (report)
398 {
399 Pout<< "mostAlignedEdge : finding nearest edge among "
400 << UIndirectList<edge>(cutFaces().edges(), pEdges)
401 << " connected to point " << pointi
402 << " coord:" << localPoints[pointi]
403 << " running between " << edgeStart << " coord:"
404 << localPoints[edgeStart]
405 << " and " << edgeEnd << " coord:"
406 << localPoints[edgeEnd]
407 << endl;
408 }
409
410 // Find the edge that gets us nearest end.
411
412 label maxEdgeI = -1;
413 scalar maxCos = -GREAT;
414
415 forAll(pEdges, i)
416 {
417 label edgeI = pEdges[i];
418
419 if
420 (
421 !(
422 patchDivision
423 && cutToMasterEdges[edgeI] == -1
424 )
425 || (
426 patchDivision
427 && regionEdge(slaveMesh, cutToSlaveEdges[edgeI])
428 )
429 )
430 {
431 const edge& e = cutFaces().edges()[edgeI];
432
433 label otherPointi = e.otherVertex(pointi);
434
435 if (otherPointi == edgeEnd)
436 {
437 // Shortcut: found edge end point.
438 if (report)
439 {
440 Pout<< " mostAlignedEdge : found end point " << edgeEnd
441 << endl;
442 }
443 return edgeI;
444 }
445
446 // Get angle between edge and edge to masterEnd
447
448 vector eVec(localPoints[otherPointi] - localPoints[pointi]);
449
450 scalar magEVec = mag(eVec);
451
452 if (magEVec < VSMALL)
453 {
455 << "Crossing zero sized edge " << edgeI
456 << " coords:" << localPoints[otherPointi]
457 << localPoints[pointi]
458 << " when walking from " << localPoints[edgeStart]
459 << " to " << localPoints[edgeEnd]
460 << endl;
461 return edgeI;
462 }
463
464 eVec /= magEVec;
465
466 const vector eToEndPoint =
468 (
469 localPoints[edgeEnd] - localPoints[otherPointi]
470 );
471
472 scalar cosAngle = eVec & eToEndPoint;
473
474 if (report)
475 {
476 Pout<< " edge:" << e << " points:" << localPoints[pointi]
477 << localPoints[otherPointi]
478 << " vec:" << eVec
479 << " vecToEnd:" << eToEndPoint
480 << " cosAngle:" << cosAngle
481 << endl;
482 }
483
484 if (cosAngle > maxCos)
485 {
486 maxCos = cosAngle;
487 maxEdgeI = edgeI;
488 }
489 }
490 }
491
492 if (maxCos > 1 - angleTol_)
493 {
494 return maxEdgeI;
495 }
496 else
497 {
498 return -1;
499 }
500}
501
502
503void Foam::faceCoupleInfo::setCutEdgeToPoints(const labelList& cutToMasterEdges)
504{
505 labelListList masterToCutEdges
506 (
508 (
509 masterPatch().nEdges(),
510 cutToMasterEdges
511 )
512 );
513
514 const edgeList& cutEdges = cutFaces().edges();
515
516 // Size extra big so searching is faster
517 cutEdgeToPoints_.reserve
518 (
519 masterPatch().nEdges()
520 + slavePatch().nEdges()
521 + cutEdges.size()
522 );
523
524 forAll(masterToCutEdges, masterEdgeI)
525 {
526 const edge& masterE = masterPatch().edges()[masterEdgeI];
527
528 //Pout<< "Master:" << masterPatch().localPoints()[masterE[0]] << ' '
529 // << masterPatch().localPoints()[masterE[1]] << endl;
530
531 const labelList& stringedEdges = masterToCutEdges[masterEdgeI];
532
533 if (stringedEdges.empty())
534 {
536 << "Did not match all of master edges to cutFace edges"
537 << nl
538 << "First unmatched edge:" << masterEdgeI << " endPoints:"
539 << masterPatch().localPoints()[masterE[0]]
540 << masterPatch().localPoints()[masterE[1]]
541 << endl
542 << "This usually means that the slave patch is not a"
543 << " subdivision of the master patch"
544 << abort(FatalError);
545 }
546 else if (stringedEdges.size() > 1)
547 {
548 // String up the edges between e[0] and e[1]. Store the points
549 // inbetween e[0] and e[1] (all in cutFaces() labels)
550
551 DynamicList<label> splitPoints(stringedEdges.size()-1);
552
553 // Unsplit edge endpoints
554 const edge unsplitEdge
555 (
556 masterToCutPoints_[masterE[0]],
557 masterToCutPoints_[masterE[1]]
558 );
559
560 label startVertI = unsplitEdge[0];
561 label startEdgeI = -1;
562
563 while (startVertI != unsplitEdge[1])
564 {
565 // Loop over all string of edges. Update
566 // - startVertI : previous vertex
567 // - startEdgeI : previous edge
568 // and insert any points into splitPoints
569
570 // For checking
571 label oldStart = startVertI;
572
573 forAll(stringedEdges, i)
574 {
575 label edgeI = stringedEdges[i];
576
577 if (edgeI != startEdgeI)
578 {
579 const edge& e = cutEdges[edgeI];
580
581 //Pout<< " cut:" << e << " at:"
582 // << cutFaces().localPoints()[e[0]]
583 // << ' ' << cutFaces().localPoints()[e[1]] << endl;
584
585 if (e[0] == startVertI)
586 {
587 startEdgeI = edgeI;
588 startVertI = e[1];
589 if (e[1] != unsplitEdge[1])
590 {
591 splitPoints.append(e[1]);
592 }
593 break;
594 }
595 else if (e[1] == startVertI)
596 {
597 startEdgeI = edgeI;
598 startVertI = e[0];
599 if (e[0] != unsplitEdge[1])
600 {
601 splitPoints.append(e[0]);
602 }
603 break;
604 }
605 }
606 }
607
608 // Check
609 if (oldStart == startVertI)
610 {
612 << " unsplitEdge:" << unsplitEdge
613 << " does not correspond to split edges "
614 << UIndirectList<edge>(cutEdges, stringedEdges)
615 << abort(FatalError);
616 }
617 }
618
619 //Pout<< "For master edge:"
620 // << unsplitEdge
621 // << " Found stringed points "
622 // << UIndirectList<point>
623 // (
624 // cutFaces().localPoints(),
625 // splitPoints.shrink()
626 // )
627 // << endl;
628
629 cutEdgeToPoints_.insert(unsplitEdge, splitPoints.shrink());
630 }
631 }
632}
633
634
635Foam::label Foam::faceCoupleInfo::matchFaces
636(
637 const scalar absTol,
638 const pointField& points0,
639 const face& f0,
640 const pointField& points1,
641 const face& f1,
642 const bool sameOrientation
643)
644{
645 if (f0.size() != f1.size())
646 {
648 << "Different sizes for supposedly matching faces." << nl
649 << "f0:" << f0 << " coords:" << UIndirectList<point>(points0, f0)
650 << nl
651 << "f1:" << f1 << " coords:" << UIndirectList<point>(points1, f1)
652 << abort(FatalError);
653 }
654
655 const scalar absTolSqr = sqr(absTol);
656
657
658 label matchFp = -1;
659
660 forAll(f0, startFp)
661 {
662 // See -if starting from startFp on f0- the two faces match.
663 bool fullMatch = true;
664
665 label fp0 = startFp;
666 label fp1 = 0;
667
668 forAll(f1, i)
669 {
670 scalar distSqr = Foam::magSqr(points0[f0[fp0]] - points1[f1[fp1]]);
671
672 if (distSqr > absTolSqr)
673 {
674 fullMatch = false;
675 break;
676 }
677
678 fp0 = f0.fcIndex(fp0); // walk forward
679
680 if (sameOrientation)
681 {
682 fp1 = f1.fcIndex(fp1);
683 }
684 else
685 {
686 fp1 = f1.rcIndex(fp1);
687 }
688 }
689
690 if (fullMatch)
691 {
692 matchFp = startFp;
693 break;
694 }
695 }
696
697 if (matchFp == -1)
698 {
700 << "No unique match between two faces" << nl
701 << "Face " << f0 << " coords "
703 << "Face " << f1 << " coords "
704 << UIndirectList<point>(points1, f1)
705 << "when using tolerance " << absTol
706 << " and forwardMatching:" << sameOrientation
707 << abort(FatalError);
708 }
709
710 return matchFp;
711}
712
713
714bool Foam::faceCoupleInfo::matchPointsThroughFaces
715(
716 const scalar absTol,
717 const pointField& cutPoints,
718 const faceList& cutFaces,
719 const pointField& patchPoints,
720 const faceList& patchFaces,
721 const bool sameOrientation,
722
723 labelList& patchToCutPoints, // patch to (uncompacted) cut points
724 labelList& cutToCompact, // compaction list for cut points
725 labelList& compactToCut // inverse ,,
726)
727{
728 // Find correspondence from patch points to cut points. This might detect
729 // shared points so the output is a patch-to-cut point list and a compaction
730 // list for the cut points (which will always be equal or more connected
731 // than the patch). Returns true if there are any duplicates.
732
733 // From slave to cut point
734 patchToCutPoints.setSize(patchPoints.size());
735 patchToCutPoints = -1;
736
737 // Compaction list for cut points: either -1 or index into master which
738 // gives the point to compact to.
739 labelList cutPointRegion(cutPoints.size(), -1);
740 DynamicList<label> cutPointRegionMaster;
741
742 forAll(patchFaces, patchFacei)
743 {
744 const face& patchF = patchFaces[patchFacei];
745
746 //const face& cutF = cutFaces[patchToCutFaces[patchFacei]];
747 const face& cutF = cutFaces[patchFacei];
748
749 // Do geometric matching to get position of cutF[0] in patchF
750 label patchFp = matchFaces
751 (
752 absTol,
753 patchPoints,
754 patchF,
755 cutPoints,
756 cutF,
757 sameOrientation // orientation
758 );
759
760 forAll(cutF, cutFp)
761 {
762 label cutPointi = cutF[cutFp];
763 label patchPointi = patchF[patchFp];
764
765 //const point& cutPt = cutPoints[cutPointi];
766 //const point& patchPt = patchPoints[patchPointi];
767 //if (mag(cutPt - patchPt) > SMALL)
768 //{
769 // FatalErrorInFunction
770 // << "cutP:" << cutPt
771 // << " patchP:" << patchPt
772 // << abort(FatalError);
773 //}
774
775 if (patchToCutPoints[patchPointi] == -1)
776 {
777 patchToCutPoints[patchPointi] = cutPointi;
778 }
779 else if (patchToCutPoints[patchPointi] != cutPointi)
780 {
781 // Multiple cut points connecting to same patch.
782 // Check if already have region & region master for this set
783 label otherCutPointi = patchToCutPoints[patchPointi];
784
785 //Pout<< "PatchPoint:" << patchPt
786 // << " matches to:" << cutPointi
787 // << " coord:" << cutPoints[cutPointi]
788 // << " and to:" << otherCutPointi
789 // << " coord:" << cutPoints[otherCutPointi]
790 // << endl;
791
792 if (cutPointRegion[otherCutPointi] != -1)
793 {
794 // Have region for this set. Copy.
795 label region = cutPointRegion[otherCutPointi];
796 cutPointRegion[cutPointi] = region;
797
798 // Update region master with min point label
799 cutPointRegionMaster[region] = min
800 (
801 cutPointRegionMaster[region],
802 cutPointi
803 );
804 }
805 else
806 {
807 // Create new region.
808 label region = cutPointRegionMaster.size();
809 cutPointRegionMaster.append
810 (
811 min(cutPointi, otherCutPointi)
812 );
813 cutPointRegion[cutPointi] = region;
814 cutPointRegion[otherCutPointi] = region;
815 }
816 }
817
818 if (sameOrientation)
819 {
820 patchFp = patchF.fcIndex(patchFp);
821 }
822 else
823 {
824 patchFp = patchF.rcIndex(patchFp);
825 }
826 }
827 }
828
829 // Rework region&master into compaction array
830 compactToCut.setSize(cutPointRegion.size());
831 cutToCompact.setSize(cutPointRegion.size());
832 cutToCompact = -1;
833 label compactPointi = 0;
834
835 forAll(cutPointRegion, i)
836 {
837 if (cutPointRegion[i] == -1)
838 {
839 // Unduplicated point. Allocate new compacted point.
840 cutToCompact[i] = compactPointi;
841 compactToCut[compactPointi] = i;
842 compactPointi++;
843 }
844 else
845 {
846 // Duplicate point. Get master.
847
848 label masterPointi = cutPointRegionMaster[cutPointRegion[i]];
849
850 if (cutToCompact[masterPointi] == -1)
851 {
852 cutToCompact[masterPointi] = compactPointi;
853 compactToCut[compactPointi] = masterPointi;
854 compactPointi++;
855 }
856 cutToCompact[i] = cutToCompact[masterPointi];
857 }
858 }
859 compactToCut.setSize(compactPointi);
860
861 return compactToCut.size() != cutToCompact.size();
862}
863
864
865Foam::scalar Foam::faceCoupleInfo::maxDistance
866(
867 const face& cutF,
868 const pointField& cutPoints,
869 const face& masterF,
870 const pointField& masterPoints
871)
872{
873 scalar maxDist = -GREAT;
874
875 forAll(cutF, fp)
876 {
877 const point& cutPt = cutPoints[cutF[fp]];
878
879 pointHit pHit = masterF.nearestPoint(cutPt, masterPoints);
880
881 maxDist = max(maxDist, pHit.distance());
882 }
883 return maxDist;
884}
885
886
887void Foam::faceCoupleInfo::findPerfectMatchingFaces
888(
889 const primitiveMesh& mesh0,
890 const primitiveMesh& mesh1,
891 const scalar absTol,
892
893 labelList& mesh0Faces,
894 labelList& mesh1Faces
895)
896{
897 // Quick check: skip face matching if either mesh has no faces
898 if (!mesh0.nFaces() || !mesh1.nFaces())
899 {
900 mesh0Faces.clear();
901 mesh1Faces.clear();
902
903 return;
904 }
905
906 // Face centres of external faces (without invoking
907 // mesh.faceCentres since mesh might have been clearedOut)
908
909 pointField fc0
910 (
911 calcFaceCentres<List>
912 (
913 mesh0.faces(),
914 mesh0.points(),
915 mesh0.nInternalFaces(),
916 mesh0.nBoundaryFaces()
917 )
918 );
919
920 pointField fc1
921 (
922 calcFaceCentres<List>
923 (
924 mesh1.faces(),
925 mesh1.points(),
926 mesh1.nInternalFaces(),
927 mesh1.nBoundaryFaces()
928 )
929 );
930
931
932 if (debug)
933 {
934 Pout<< "Face matching tolerance : " << absTol << endl;
935 }
936
937
938 // Match geometrically
939 labelList from1To0;
940 bool matchedAllFaces = matchPoints
941 (
942 fc1,
943 fc0,
944 scalarField(fc1.size(), absTol),
945 false,
946 from1To0
947 );
948
949 if (matchedAllFaces)
950 {
952 << "Matched ALL " << fc1.size()
953 << " boundary faces of mesh0 to boundary faces of mesh1." << endl
954 << "This is only valid if the mesh to add is fully"
955 << " enclosed by the mesh it is added to." << endl;
956 }
957
958
959 // Collect matches.
960 label nMatched = 0;
961
962 mesh0Faces.setSize(fc0.size());
963 mesh1Faces.setSize(fc1.size());
964
965 forAll(from1To0, i)
966 {
967 if (from1To0[i] != -1)
968 {
969 mesh1Faces[nMatched] = i + mesh1.nInternalFaces();
970 mesh0Faces[nMatched] = from1To0[i] + mesh0.nInternalFaces();
971
972 nMatched++;
973 }
974 }
975
976 mesh0Faces.setSize(nMatched);
977 mesh1Faces.setSize(nMatched);
978}
979
980
981void Foam::faceCoupleInfo::findSlavesCoveringMaster
982(
983 const primitiveMesh& mesh0,
984 const primitiveMesh& mesh1,
985 const scalar absTol,
986
987 labelList& mesh0Faces,
988 labelList& mesh1Faces
989)
990{
991 // Construct octree from all mesh0 boundary faces
992
993 Random rndGen(123456);
994
995 treeBoundBox overallBb(mesh0.points());
996
998 (
1000 (
1001 mesh0,
1002 // boundary faces only
1003 identity(mesh0.nBoundaryFaces(), mesh0.nInternalFaces())
1004 ),
1005
1006 overallBb.extend(rndGen, 1e-4), // overall search domain
1007 8, // maxLevel
1008 10, // leafsize
1009 3.0 // duplicity
1010 );
1011
1012 if (debug)
1013 {
1014 Pout<< "findSlavesCoveringMaster :"
1015 << " constructed octree for mesh0 boundary faces" << endl;
1016 }
1017
1018 // Search nearest face for every mesh1 boundary face.
1019
1020 labelHashSet mesh0Set(mesh0.nBoundaryFaces());
1021 labelHashSet mesh1Set(mesh1.nBoundaryFaces());
1022
1023 for
1024 (
1025 label mesh1Facei = mesh1.nInternalFaces();
1026 mesh1Facei < mesh1.nFaces();
1027 mesh1Facei++
1028 )
1029 {
1030 const face& f1 = mesh1.faces()[mesh1Facei];
1031
1032 // Generate face centre (prevent cellCentres() reconstruction)
1033 point fc(f1.centre(mesh1.points()));
1034
1035 pointIndexHit nearInfo = tree.findNearest(fc, Foam::sqr(absTol));
1036
1037 if (nearInfo.hit())
1038 {
1039 label mesh0Facei = tree.shapes().objectIndex(nearInfo.index());
1040
1041 // Check if points of f1 actually lie on top of mesh0 face
1042 // This is the bit that might fail since if f0 is severely warped
1043 // and f1's points are not present on f0 (since subdivision)
1044 // then the distance of the points to f0 might be quite large
1045 // and the test will fail. This all should in fact be some kind
1046 // of walk from known corresponding points/edges/faces.
1047 scalar dist =
1048 maxDistance
1049 (
1050 f1,
1051 mesh1.points(),
1052 mesh0.faces()[mesh0Facei],
1053 mesh0.points()
1054 );
1055
1056 if (dist < absTol)
1057 {
1058 mesh0Set.insert(mesh0Facei);
1059 mesh1Set.insert(mesh1Facei);
1060 }
1061 }
1062 }
1063
1064 if (debug)
1065 {
1066 Pout<< "findSlavesCoveringMaster :"
1067 << " matched " << mesh1Set.size() << " mesh1 faces to "
1068 << mesh0Set.size() << " mesh0 faces" << endl;
1069 }
1070
1071 mesh0Faces = mesh0Set.toc();
1072 mesh1Faces = mesh1Set.toc();
1073}
1074
1075
1076Foam::label Foam::faceCoupleInfo::growCutFaces
1077(
1078 const labelList& cutToMasterEdges,
1079 Map<labelList>& candidates
1080)
1081{
1082 if (debug)
1083 {
1084 Pout<< "growCutFaces :"
1085 << " growing cut faces to masterPatch" << endl;
1086 }
1087
1088 label nTotChanged = 0;
1089
1090 while (true)
1091 {
1092 const labelListList& cutFaceEdges = cutFaces().faceEdges();
1093 const labelListList& cutEdgeFaces = cutFaces().edgeFaces();
1094
1095 label nChanged = 0;
1096
1097 forAll(cutToMasterFaces_, cutFacei)
1098 {
1099 const label masterFacei = cutToMasterFaces_[cutFacei];
1100
1101 if (masterFacei != -1)
1102 {
1103 // We now have a cutFace for which we have already found a
1104 // master face. Grow this masterface across any internal edge
1105 // (internal: no corresponding master edge)
1106
1107 const labelList& fEdges = cutFaceEdges[cutFacei];
1108
1109 forAll(fEdges, i)
1110 {
1111 const label cutEdgeI = fEdges[i];
1112
1113 if (cutToMasterEdges[cutEdgeI] == -1)
1114 {
1115 // So this edge:
1116 // - internal to the cutPatch (since all region edges
1117 // marked before)
1118 // - on cutFacei which corresponds to masterFace.
1119 // Mark all connected faces with this masterFace.
1120
1121 const labelList& eFaces = cutEdgeFaces[cutEdgeI];
1122
1123 forAll(eFaces, j)
1124 {
1125 const label facei = eFaces[j];
1126
1127 if (cutToMasterFaces_[facei] == -1)
1128 {
1129 cutToMasterFaces_[facei] = masterFacei;
1130 candidates.erase(facei);
1131 nChanged++;
1132 }
1133 else if (cutToMasterFaces_[facei] != masterFacei)
1134 {
1135 const pointField& cutPoints =
1136 cutFaces().points();
1137 const pointField& masterPoints =
1138 masterPatch().points();
1139
1140 const edge& e = cutFaces().edges()[cutEdgeI];
1141
1142 label myMaster = cutToMasterFaces_[facei];
1143 const face& myF = masterPatch()[myMaster];
1144
1145 const face& nbrF = masterPatch()[masterFacei];
1146
1148 << "Cut face "
1149 << cutFaces()[facei].points(cutPoints)
1150 << " has master "
1151 << myMaster
1152 << " but also connects to nbr face "
1153 << cutFaces()[cutFacei].points(cutPoints)
1154 << " with master " << masterFacei
1155 << nl
1156 << "myMasterFace:"
1157 << myF.points(masterPoints)
1158 << " nbrMasterFace:"
1159 << nbrF.points(masterPoints) << nl
1160 << "Across cut edge " << cutEdgeI
1161 << " coord:"
1162 << cutFaces().localPoints()[e[0]]
1163 << cutFaces().localPoints()[e[1]]
1164 << abort(FatalError);
1165 }
1166 }
1167 }
1168 }
1169 }
1170 }
1171
1172 if (debug)
1173 {
1174 Pout<< "growCutFaces : Grown an additional " << nChanged
1175 << " cut to master face" << " correspondence" << endl;
1176 }
1177
1178 nTotChanged += nChanged;
1179
1180 if (nChanged == 0)
1181 {
1182 break;
1183 }
1184 }
1185
1186 return nTotChanged;
1187}
1188
1189
1190void Foam::faceCoupleInfo::checkMatch(const labelList& cutToMasterEdges) const
1191{
1192 const pointField& cutLocalPoints = cutFaces().localPoints();
1193
1194 const pointField& masterLocalPoints = masterPatch().localPoints();
1195 const faceList& masterLocalFaces = masterPatch().localFaces();
1196
1197 forAll(cutToMasterEdges, cutEdgeI)
1198 {
1199 const edge& e = cutFaces().edges()[cutEdgeI];
1200
1201 if (cutToMasterEdges[cutEdgeI] == -1)
1202 {
1203 // Internal edge. Check that master face is same on both sides.
1204 const labelList& cutEFaces = cutFaces().edgeFaces()[cutEdgeI];
1205
1206 label masterFacei = -1;
1207
1208 forAll(cutEFaces, i)
1209 {
1210 label cutFacei = cutEFaces[i];
1211
1212 if (cutToMasterFaces_[cutFacei] != -1)
1213 {
1214 if (masterFacei == -1)
1215 {
1216 masterFacei = cutToMasterFaces_[cutFacei];
1217 }
1218 else if (masterFacei != cutToMasterFaces_[cutFacei])
1219 {
1220 label myMaster = cutToMasterFaces_[cutFacei];
1221 const face& myF = masterLocalFaces[myMaster];
1222
1223 const face& nbrF = masterLocalFaces[masterFacei];
1224
1226 << "Internal CutEdge " << e
1227 << " coord:"
1228 << cutLocalPoints[e[0]]
1229 << cutLocalPoints[e[1]]
1230 << " connects to master " << myMaster
1231 << " and to master " << masterFacei << nl
1232 << "myMasterFace:"
1233 << myF.points(masterLocalPoints)
1234 << " nbrMasterFace:"
1235 << nbrF.points(masterLocalPoints)
1236 << abort(FatalError);
1237 }
1238 }
1239 }
1240 }
1241 }
1242}
1243
1244
1245Foam::label Foam::faceCoupleInfo::matchEdgeFaces
1246(
1247 const labelList& cutToMasterEdges,
1248 Map<labelList>& candidates
1249)
1250{
1251 // Extends matching information by elimination across cutFaces using more
1252 // than one region edge. Updates cutToMasterFaces_ and sets candidates which
1253 // is for every cutface on a region edge the possible master faces.
1254
1255 // For every unassigned cutFacei the possible list of master faces.
1256 candidates.clear();
1257 candidates.reserve(cutFaces().size());
1258
1259 label nChanged = 0;
1260
1261 forAll(cutToMasterEdges, cutEdgeI)
1262 {
1263 label masterEdgeI = cutToMasterEdges[cutEdgeI];
1264
1265 if (masterEdgeI != -1)
1266 {
1267 // So cutEdgeI is matched to masterEdgeI. For all cut faces using
1268 // the cut edge store the master face as a candidate match.
1269
1270 const labelList& cutEFaces = cutFaces().edgeFaces()[cutEdgeI];
1271 const labelList& masterEFaces =
1272 masterPatch().edgeFaces()[masterEdgeI];
1273
1274 forAll(cutEFaces, i)
1275 {
1276 label cutFacei = cutEFaces[i];
1277
1278 if (cutToMasterFaces_[cutFacei] == -1)
1279 {
1280 // So this face (cutFacei) is on an edge (cutEdgeI) that
1281 // is used by some master faces (masterEFaces). Check if
1282 // the cutFacei and some of these masterEFaces share more
1283 // than one edge (which uniquely defines face)
1284
1285 // Combine master faces with current set of candidate
1286 // master faces.
1287 auto fnd = candidates.find(cutFacei);
1288
1289 if (!fnd.good())
1290 {
1291 // No info yet for cutFacei. Add all master faces as
1292 // candidates
1293 candidates.insert(cutFacei, masterEFaces);
1294 }
1295 else
1296 {
1297 // From some other cutEdgeI there are already some
1298 // candidate master faces. Check the overlap with
1299 // the current set of master faces.
1300 const labelList& masterFaces = fnd.val();
1301
1302 DynamicList<label> newCandidates(masterFaces.size());
1303
1304 forAll(masterEFaces, j)
1305 {
1306 if (masterFaces.found(masterEFaces[j]))
1307 {
1308 newCandidates.append(masterEFaces[j]);
1309 }
1310 }
1311
1312 if (newCandidates.size() == 1)
1313 {
1314 // Perfect match. Delete entry from candidates map.
1315 cutToMasterFaces_[cutFacei] = newCandidates[0];
1316 candidates.erase(cutFacei);
1317 nChanged++;
1318 }
1319 else
1320 {
1321 // Should not happen?
1322 //Pout<< "On edge:" << cutEdgeI
1323 // << " have connected masterFaces:"
1324 // << masterEFaces
1325 // << " and from previous edge we have"
1326 // << " connected masterFaces" << masterFaces
1327 // << " . Overlap:" << newCandidates << endl;
1328
1329 fnd() = newCandidates.shrink();
1330 }
1331 }
1332 }
1333
1334 }
1335 }
1336 }
1337
1338 if (debug)
1339 {
1340 Pout<< "matchEdgeFaces : Found " << nChanged
1341 << " faces where there was"
1342 << " only one remaining choice for cut-master correspondence"
1343 << endl;
1344 }
1345
1346 return nChanged;
1347}
1348
1349
1350Foam::label Foam::faceCoupleInfo::geometricMatchEdgeFaces
1351(
1352 Map<labelList>& candidates
1353)
1354{
1355 const pointField& cutPoints = cutFaces().points();
1356
1357 label nChanged = 0;
1358
1359 // Mark all master faces that have been matched so far.
1360
1361 labelListList masterToCutFaces
1362 (
1364 (
1365 masterPatch().size(),
1366 cutToMasterFaces_
1367 )
1368 );
1369
1370 forAllConstIters(candidates, iter)
1371 {
1372 label cutFacei = iter.key();
1373 const labelList& masterFaces = iter.val();
1374
1375 const face& cutF = cutFaces()[cutFacei];
1376
1377 if (cutToMasterFaces_[cutFacei] == -1)
1378 {
1379 // Find the best matching master face.
1380 scalar minDist = GREAT;
1381 label minMasterFacei = -1;
1382
1383 forAll(masterFaces, i)
1384 {
1385 label masterFacei = masterFaces[i];
1386
1387 if (masterToCutFaces[masterFaces[i]].empty())
1388 {
1389 scalar dist = maxDistance
1390 (
1391 cutF,
1392 cutPoints,
1393 masterPatch()[masterFacei],
1394 masterPatch().points()
1395 );
1396
1397 if (dist < minDist)
1398 {
1399 minDist = dist;
1400 minMasterFacei = masterFacei;
1401 }
1402 }
1403 }
1404
1405 if (minMasterFacei != -1)
1406 {
1407 cutToMasterFaces_[cutFacei] = minMasterFacei;
1408 masterToCutFaces[minMasterFacei] = cutFacei;
1409 nChanged++;
1410 }
1411 }
1412 }
1413
1414 // (inefficiently) force candidates to be uptodate.
1415 forAll(cutToMasterFaces_, cutFacei)
1416 {
1417 if (cutToMasterFaces_[cutFacei] != -1)
1418 {
1419 candidates.erase(cutFacei);
1420 }
1421 }
1422
1423
1424 if (debug)
1425 {
1426 Pout<< "geometricMatchEdgeFaces : Found " << nChanged
1427 << " faces where there was"
1428 << " only one remaining choice for cut-master correspondence"
1429 << endl;
1430 }
1431
1432 return nChanged;
1433}
1434
1435
1436void Foam::faceCoupleInfo::perfectPointMatch
1437(
1438 const scalar absTol,
1439 const bool slaveFacesOrdered
1440)
1441{
1442 // Calculate the set of cut faces inbetween master and slave patch assuming
1443 // perfect match (and optional face ordering on slave)
1444
1445 if (debug)
1446 {
1447 Pout<< "perfectPointMatch :"
1448 << " Matching master and slave to cut."
1449 << " Master and slave faces are identical" << nl;
1450
1451 if (slaveFacesOrdered)
1452 {
1453 Pout<< "and master and slave faces are ordered"
1454 << " (on coupled patches)" << endl;
1455 }
1456 else
1457 {
1458 Pout<< "and master and slave faces are not ordered"
1459 << " (on coupled patches)" << endl;
1460 }
1461 }
1462
1463 cutToMasterFaces_ = identity(masterPatch().size());
1464 cutPoints_ = masterPatch().localPoints();
1465 cutFacesPtr_.reset
1466 (
1468 (
1469 masterPatch().localFaces(),
1470 cutPoints_
1471 )
1472 );
1473 masterToCutPoints_ = identity(cutPoints_.size());
1474
1475
1476 // Cut faces to slave patch.
1477 bool matchedAllFaces = false;
1478
1479 if (slaveFacesOrdered)
1480 {
1481 cutToSlaveFaces_ = identity(cutFaces().size());
1482 matchedAllFaces = (cutFaces().size() == slavePatch().size());
1483 }
1484 else
1485 {
1486 // Faces do not have to be ordered (but all have
1487 // to match). Note: Faces will be already ordered if we enter here from
1488 // construct from meshes.
1489
1490 matchedAllFaces = matchPoints
1491 (
1492 calcFaceCentres<List>
1493 (
1494 cutFaces(),
1495 cutFaces().points(),
1496 0,
1497 cutFaces().size()
1498 ),
1499 calcFaceCentres<IndirectList>
1500 (
1501 slavePatch(),
1502 slavePatch().points(),
1503 0,
1504 slavePatch().size()
1505 ),
1506 scalarField(slavePatch().size(), absTol),
1507 false,
1508 cutToSlaveFaces_
1509 );
1510
1511 // If some of the face centres did not match, then try to match the
1512 // point averages instead. There is no division by the face area in
1513 // calculating the point average, so this is more stable when faces
1514 // collapse onto a line or point.
1515 if (!matchedAllFaces)
1516 {
1517 labelList cutToSlaveFacesTemp(cutToSlaveFaces_.size(), -1);
1518
1520 (
1521 calcFacePointAverages<List>
1522 (
1523 cutFaces(),
1524 cutFaces().points(),
1525 0,
1526 cutFaces().size()
1527 ),
1528 calcFacePointAverages<IndirectList>
1529 (
1530 slavePatch(),
1531 slavePatch().points(),
1532 0,
1533 slavePatch().size()
1534 ),
1535 scalarField(slavePatch().size(), absTol),
1536 true,
1537 cutToSlaveFacesTemp
1538 );
1539
1540 cutToSlaveFaces_ = max(cutToSlaveFaces_, cutToSlaveFacesTemp);
1541
1542 matchedAllFaces = min(cutToSlaveFaces_) != -1;
1543 }
1544 }
1545
1546
1547 if (!matchedAllFaces)
1548 {
1550 << "Did not match all of the master faces to the slave faces"
1551 << endl
1552 << "This usually means that the slave patch and master patch"
1553 << " do not align to within " << absTol << " metre."
1554 << abort(FatalError);
1555 }
1556
1557
1558 // Find correspondence from slave points to cut points. This might
1559 // detect shared points so the output is a slave-to-cut point list
1560 // and a compaction list.
1561
1562 labelList cutToCompact, compactToCut;
1563 matchPointsThroughFaces
1564 (
1565 absTol,
1566 cutFaces().localPoints(),
1567 reorder(cutToSlaveFaces_, cutFaces().localFaces()),
1568 slavePatch().localPoints(),
1569 slavePatch().localFaces(),
1570 false, // slave and cut have opposite orientation
1571
1572 slaveToCutPoints_, // slave to (uncompacted) cut points
1573 cutToCompact, // compaction map: from cut to compacted
1574 compactToCut // compaction map: from compacted to cut
1575 );
1576
1577
1578 // Use compaction lists to renumber cutPoints.
1579 cutPoints_ = UIndirectList<point>(cutPoints_, compactToCut)();
1580 {
1581 const faceList& cutLocalFaces = cutFaces().localFaces();
1582
1583 faceList compactFaces(cutLocalFaces.size());
1584 forAll(cutLocalFaces, i)
1585 {
1586 compactFaces[i] = renumber(cutToCompact, cutLocalFaces[i]);
1587 }
1588 cutFacesPtr_.reset
1589 (
1591 (
1592 compactFaces,
1593 cutPoints_
1594 )
1595 );
1596 }
1597 inplaceRenumber(cutToCompact, slaveToCutPoints_);
1598 inplaceRenumber(cutToCompact, masterToCutPoints_);
1599}
1600
1601
1602void Foam::faceCoupleInfo::subDivisionMatch
1603(
1604 const polyMesh& slaveMesh,
1605 const bool patchDivision,
1606 const scalar absTol
1607)
1608{
1609 if (debug)
1610 {
1611 Pout<< "subDivisionMatch :"
1612 << " Matching master and slave to cut."
1613 << " Slave can be subdivision of master but all master edges"
1614 << " have to be covered by slave edges." << endl;
1615 }
1616
1617 // Assume slave patch is subdivision of the master patch so cutFaces is a
1618 // copy of the slave (but reversed since orientation has to be like master).
1619 cutPoints_ = slavePatch().localPoints();
1620 {
1621 faceList cutFaces(slavePatch().size());
1622
1623 forAll(cutFaces, i)
1624 {
1625 cutFaces[i] = slavePatch().localFaces()[i].reverseFace();
1626 }
1627 cutFacesPtr_.reset(new primitiveFacePatch(cutFaces, cutPoints_));
1628 }
1629
1630 // Cut is copy of slave so addressing to slave is trivial.
1631 cutToSlaveFaces_ = identity(cutFaces().size());
1632 slaveToCutPoints_ = identity(slavePatch().nPoints());
1633
1634 // Cut edges to slave patch
1635 labelList cutToSlaveEdges
1636 (
1637 findMappedEdges
1638 (
1639 cutFaces().edges(),
1640 slaveToCutPoints_, //note:should be cutToSlavePoints but since iden
1641 slavePatch()
1642 )
1643 );
1644
1645
1646 if (debug)
1647 {
1648 OFstream str("regionEdges.obj");
1649
1650 Pout<< "subDivisionMatch :"
1651 << " addressing for slave patch fully done."
1652 << " Dumping region edges to " << str.name() << endl;
1653
1654 label vertI = 0;
1655
1656 forAll(slavePatch().edges(), slaveEdgeI)
1657 {
1658 if (regionEdge(slaveMesh, slaveEdgeI))
1659 {
1660 const edge& e = slavePatch().edges()[slaveEdgeI];
1661
1662 meshTools::writeOBJ(str, slavePatch().localPoints()[e[0]]);
1663 vertI++;
1664 meshTools::writeOBJ(str, slavePatch().localPoints()[e[1]]);
1665 vertI++;
1666 str<< "l " << vertI-1 << ' ' << vertI << nl;
1667 }
1668 }
1669 }
1670
1671
1672 // Addressing from cut to master.
1673
1674 // Cut points to master patch
1675 // - determine master points to cut points. Has to be full!
1676 // - invert to get cut to master
1677 if (debug)
1678 {
1679 Pout<< "subDivisionMatch :"
1680 << " matching master points to cut points" << endl;
1681 }
1682
1683 bool matchedAllPoints = matchPoints
1684 (
1685 masterPatch().localPoints(),
1686 cutPoints_,
1687 scalarField(masterPatch().nPoints(), absTol),
1688 true,
1689 masterToCutPoints_
1690 );
1691
1692 if (!matchedAllPoints)
1693 {
1695 << "Did not match all of the master points to the slave points"
1696 << endl
1697 << "This usually means that the slave patch is not a"
1698 << " subdivision of the master patch"
1699 << abort(FatalError);
1700 }
1701
1702
1703 // Do masterEdges to cutEdges :
1704 // - mark all edges between two masterEdge endpoints. (geometric test since
1705 // this is the only distinction)
1706 // - this gives the borders inbetween which all cutfaces come from
1707 // a single master face.
1708 if (debug)
1709 {
1710 Pout<< "subDivisionMatch :"
1711 << " matching cut edges to master edges" << endl;
1712 }
1713
1714 const edgeList& masterEdges = masterPatch().edges();
1715 const pointField& masterPoints = masterPatch().localPoints();
1716
1717 const edgeList& cutEdges = cutFaces().edges();
1718
1719 labelList cutToMasterEdges(cutFaces().nEdges(), -1);
1720
1721 forAll(masterEdges, masterEdgeI)
1722 {
1723 const edge& masterEdge = masterEdges[masterEdgeI];
1724
1725 label cutPoint0 = masterToCutPoints_[masterEdge[0]];
1726 label cutPoint1 = masterToCutPoints_[masterEdge[1]];
1727
1728 // Find edges between cutPoint0 and cutPoint1.
1729
1730 label cutPointi = cutPoint0;
1731 do
1732 {
1733 // Find edge (starting at pointi on cut), aligned with master
1734 // edge.
1735 label cutEdgeI =
1736 mostAlignedCutEdge
1737 (
1738 false,
1739 slaveMesh,
1740 patchDivision,
1741 cutToMasterEdges,
1742 cutToSlaveEdges,
1743 cutPointi,
1744 cutPoint0,
1745 cutPoint1
1746 );
1747
1748 if (cutEdgeI == -1)
1749 {
1750 // Problem. Report while matching to produce nice error message
1751 mostAlignedCutEdge
1752 (
1753 true,
1754 slaveMesh,
1755 patchDivision,
1756 cutToMasterEdges,
1757 cutToSlaveEdges,
1758 cutPointi,
1759 cutPoint0,
1760 cutPoint1
1761 );
1762
1763 Pout<< "Dumping unmatched pointEdges to errorEdges.obj"
1764 << endl;
1765
1766 writeOBJ
1767 (
1768 "errorEdges.obj",
1769 edgeList
1770 (
1772 (
1773 cutFaces().edges(),
1774 cutFaces().pointEdges()[cutPointi]
1775 )
1776 ),
1777 cutFaces().localPoints(),
1778 false
1779 );
1780
1782 << "Problem in finding cut edges corresponding to"
1783 << " master edge " << masterEdge
1784 << " points:" << masterPoints[masterEdge[0]]
1785 << ' ' << masterPoints[masterEdge[1]]
1786 << " corresponding cut edge: (" << cutPoint0
1787 << ") " << cutPoint1
1788 << abort(FatalError);
1789 }
1790
1791 cutToMasterEdges[cutEdgeI] = masterEdgeI;
1792
1793 cutPointi = cutEdges[cutEdgeI].otherVertex(cutPointi);
1794
1795 } while (cutPointi != cutPoint1);
1796 }
1797
1798 if (debug)
1799 {
1800 // Write all matched edges.
1801 writeEdges(cutToMasterEdges, cutToSlaveEdges);
1802 }
1803
1804 // Rework cutToMasterEdges into list of points inbetween two endpoints
1805 // (cutEdgeToPoints_)
1806 setCutEdgeToPoints(cutToMasterEdges);
1807
1808
1809 // Now that we have marked the cut edges that lie on top of master edges
1810 // we can find cut faces that have two (or more) of these edges.
1811 // Note that there can be multiple faces having two or more matched edges
1812 // since the cut faces can form a non-manifold surface(?)
1813 // So the matching is done as an elimination process where for every
1814 // cutFace on a matched edge we store the possible master faces and
1815 // eliminate more and more until we only have one possible master face
1816 // left.
1817
1818 if (debug)
1819 {
1820 Pout<< "subDivisionMatch :"
1821 << " matching (topological) cut faces to masterPatch" << endl;
1822 }
1823
1824 // For every unassigned cutFacei the possible list of master faces.
1825 Map<labelList> candidates(cutFaces().size());
1826
1827 cutToMasterFaces_.setSize(cutFaces().size());
1828 cutToMasterFaces_ = -1;
1829
1830 while (true)
1831 {
1832 // See if there are any edges left where there is only one remaining
1833 // candidate.
1834 label nChanged = matchEdgeFaces(cutToMasterEdges, candidates);
1835
1836 checkMatch(cutToMasterEdges);
1837
1838
1839 // Now we should have found a face correspondence for every cutFace
1840 // that uses two or more cutEdges. Floodfill from these across
1841 // 'internal' cutedges (i.e. edges that do not have a master
1842 // equivalent) to determine some more cutToMasterFaces
1843 nChanged += growCutFaces(cutToMasterEdges, candidates);
1844
1845 checkMatch(cutToMasterEdges);
1846
1847 if (nChanged == 0)
1848 {
1849 break;
1850 }
1851 }
1852
1853 if (debug)
1854 {
1855 Pout<< "subDivisionMatch :"
1856 << " matching (geometric) cut faces to masterPatch" << endl;
1857 }
1858
1859 // Above should have matched all cases fully. Cannot prove this yet so
1860 // do any remaining unmatched edges by a geometric test.
1861 while (true)
1862 {
1863 // See if there are any edges left where there is only one remaining
1864 // candidate.
1865 label nChanged = geometricMatchEdgeFaces(candidates);
1866
1867 checkMatch(cutToMasterEdges);
1868
1869 nChanged += growCutFaces(cutToMasterEdges, candidates);
1870
1871 checkMatch(cutToMasterEdges);
1872
1873 if (nChanged == 0)
1874 {
1875 break;
1876 }
1877 }
1878
1879
1880 // All cut faces matched?
1881 forAll(cutToMasterFaces_, cutFacei)
1882 {
1883 if (cutToMasterFaces_[cutFacei] == -1)
1884 {
1885 const face& cutF = cutFaces()[cutFacei];
1886
1888 << "Did not match all of cutFaces to a master face" << nl
1889 << "First unmatched cut face:" << cutFacei << " with points:"
1890 << UIndirectList<point>(cutFaces().points(), cutF)
1891 << nl
1892 << "This usually means that the slave patch is not a"
1893 << " subdivision of the master patch"
1894 << abort(FatalError);
1895 }
1896 }
1897
1898 if (debug)
1899 {
1900 Pout<< "subDivisionMatch :"
1901 << " finished matching master and slave to cut" << endl;
1902 }
1903
1904 if (debug)
1905 {
1906 writePointsFaces();
1907 }
1908}
1909
1910
1911// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1912
1914(
1915 const polyMesh& masterMesh,
1916 const polyMesh& slaveMesh,
1917 const scalar absTol,
1918 const bool perfectMatch
1919)
1920:
1921 masterPatchPtr_(nullptr),
1922 slavePatchPtr_(nullptr),
1923 cutPoints_(0),
1924 cutFacesPtr_(nullptr),
1925 cutToMasterFaces_(0),
1926 masterToCutPoints_(0),
1927 cutToSlaveFaces_(0),
1928 slaveToCutPoints_(0),
1929 cutEdgeToPoints_(0)
1930{
1931 // Get faces on both meshes that are aligned.
1932 // (not ordered i.e. masterToMesh[0] does
1933 // not couple to slaveToMesh[0]. This ordering is done later on)
1934 labelList masterToMesh;
1935 labelList slaveToMesh;
1936
1937 if (perfectMatch)
1938 {
1939 // Identical faces so use tight face-centre comparison.
1940 findPerfectMatchingFaces
1941 (
1942 masterMesh,
1943 slaveMesh,
1944 absTol,
1945 masterToMesh,
1946 slaveToMesh
1947 );
1948 }
1949 else
1950 {
1951 // Slave subdivision of master so use 'nearest'. Bit dodgy, especially
1952 // with using absTol (which is quite small)
1953 findSlavesCoveringMaster
1954 (
1955 masterMesh,
1956 slaveMesh,
1957 absTol,
1958 masterToMesh,
1959 slaveToMesh
1960 );
1961 }
1962
1963 // Construct addressing engines for both sides
1964 masterPatchPtr_.reset
1965 (
1967 (
1968 IndirectList<face>(masterMesh.faces(), masterToMesh),
1969 masterMesh.points()
1970 )
1971 );
1972
1973 slavePatchPtr_.reset
1974 (
1976 (
1977 IndirectList<face>(slaveMesh.faces(), slaveToMesh),
1978 slaveMesh.points()
1979 )
1980 );
1981
1982
1983 if (perfectMatch)
1984 {
1985 // Faces are perfectly aligned but probably not ordered.
1986 perfectPointMatch(absTol, false);
1987 }
1988 else
1989 {
1990 // Slave faces are subdivision of master face. Faces not ordered.
1991 subDivisionMatch(slaveMesh, false, absTol);
1992 }
1993
1995 {
1996 writePointsFaces();
1997 }
1998}
1999
2000
2002(
2003 const polyMesh& masterMesh,
2004 const labelList& masterAddressing,
2005 const polyMesh& slaveMesh,
2006 const labelList& slaveAddressing,
2007 const scalar absTol,
2008 const bool perfectMatch,
2009 const bool orderedFaces,
2010 const bool patchDivision
2011)
2012:
2013 masterPatchPtr_
2014 (
2016 (
2017 IndirectList<face>(masterMesh.faces(), masterAddressing),
2018 masterMesh.points()
2019 )
2020 ),
2021 slavePatchPtr_
2022 (
2024 (
2025 IndirectList<face>(slaveMesh.faces(), slaveAddressing),
2026 slaveMesh.points()
2027 )
2028 ),
2029 cutPoints_(0),
2030 cutFacesPtr_(nullptr),
2031 cutToMasterFaces_(0),
2032 masterToCutPoints_(0),
2033 cutToSlaveFaces_(0),
2034 slaveToCutPoints_(0),
2035 cutEdgeToPoints_(0)
2036{
2037 if (perfectMatch && (masterAddressing.size() != slaveAddressing.size()))
2038 {
2039 FatalErrorInFunction
2040 << "Perfect match specified but number of master and slave faces"
2041 << " differ." << endl
2042 << "master:" << masterAddressing.size()
2043 << " slave:" << slaveAddressing.size()
2044 << abort(FatalError);
2045 }
2046
2047 if
2048 (
2049 masterAddressing.size()
2050 && min(masterAddressing) < masterMesh.nInternalFaces()
2051 )
2052 {
2053 FatalErrorInFunction
2054 << "Supplied internal face on master mesh to couple." << nl
2055 << "Faces to be coupled have to be boundary faces."
2056 << abort(FatalError);
2057 }
2058 if
2059 (
2060 slaveAddressing.size()
2061 && min(slaveAddressing) < slaveMesh.nInternalFaces()
2062 )
2063 {
2064 FatalErrorInFunction
2065 << "Supplied internal face on slave mesh to couple." << nl
2066 << "Faces to be coupled have to be boundary faces."
2067 << abort(FatalError);
2068 }
2069
2070
2071 if (perfectMatch)
2072 {
2073 perfectPointMatch(absTol, orderedFaces);
2074 }
2075 else
2076 {
2077 // Slave faces are subdivision of master face. Faces not ordered.
2078 subDivisionMatch(slaveMesh, patchDivision, absTol);
2079 }
2080
2081 if (debug)
2082 {
2083 writePointsFaces();
2084 }
2085}
2086
2087
2088// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2089
2091{
2092 labelList faces(pp.size());
2093
2094 label facei = pp.start();
2095
2096 forAll(pp, i)
2098 faces[i] = facei++;
2099 }
2100 return faces;
2101}
2102
2103
2105{
2106 Map<label> map(lst.size());
2107
2108 forAll(lst, i)
2109 {
2110 if (lst[i] != -1)
2111 {
2112 map.insert(i, lst[i]);
2113 }
2114 }
2115 return map;
2116}
2117
2118
2119Foam::Map<Foam::labelList> Foam::faceCoupleInfo::makeMap
2120(
2121 const labelListList& lst
2122)
2123{
2124 Map<labelList> map(lst.size());
2125
2126 forAll(lst, i)
2127 {
2128 if (lst[i].size())
2129 {
2130 map.insert(i, lst[i]);
2131 }
2132 }
2133 return map;
2134}
2135
2136
2137// ************************************************************************* //
if(maxValue - minValue< SMALL)
if(patchID !=-1)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
A List with indirect addressing.
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
Random number generator.
Definition Random.H:56
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
static labelList faceLabels(const polyPatch &)
Utility functions.
static Map< label > makeMap(const labelList &)
Create Map from List.
faceCoupleInfo(const polyMesh &mesh0, const polyMesh &mesh1, const scalar absTol, const bool perfectMatch)
Construct from two meshes and absolute tolerance.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
A class for handling file names.
Definition fileName.H:75
Non-pointer based hierarchical recursive searching.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Cell-face mesh analysis engine.
label nInternalFaces() const noexcept
Number of internal faces.
Standard boundBox with extra functionality for use in octree.
Encapsulation of data for searching on faces.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Determine correspondence between points. See below.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition meshTools.C:352
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
PointHit< point > pointHit
A PointHit with a 3D point.
Definition pointHit.H:51
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition matchPoints.C:29
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition ListOps.C:125
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
ListType reorder(const labelUList &oldToNew, const ListType &input, const bool prune=false)
Reorder the elements of a list.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label newPointi
Tree tree(triangles.begin(), triangles.end())
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Random rndGen
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))