FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

faceCoupleInfo.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*  \
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*----------------------------------------------------------------------------*/
00025 
00026 #include "faceCoupleInfo.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <OpenFOAM/matchPoints.H>
00029 #include <OpenFOAM/indirectPrimitivePatch.H>
00030 #include <meshTools/meshTools.H>
00031 #include <meshTools/treeDataFace.H>
00032 #include <meshTools/indexedOctree.H>
00033 #include <OpenFOAM/OFstream.H>
00034 
00035 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00036 
00037 defineTypeNameAndDebug(Foam::faceCoupleInfo, 0);
00038 
00039 const Foam::scalar Foam::faceCoupleInfo::angleTol_ = 1E-3;
00040 
00041 
00042 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00043 
00044 //- Write edges
00045 void Foam::faceCoupleInfo::writeOBJ
00046 (
00047     const fileName& fName,
00048     const edgeList& edges,
00049     const pointField& points,
00050     const bool compact
00051 )
00052 {
00053     OFstream str(fName);
00054 
00055     labelList pointMap(points.size(), -1);
00056 
00057     if (compact)
00058     {
00059         label newPointI = 0;
00060 
00061         forAll(edges, edgeI)
00062         {
00063             const edge& e = edges[edgeI];
00064 
00065             forAll(e, eI)
00066             {
00067                 label pointI = e[eI];
00068 
00069                 if (pointMap[pointI] == -1)
00070                 {
00071                     pointMap[pointI] = newPointI++;
00072 
00073                     meshTools::writeOBJ(str, points[pointI]);
00074                 }
00075             }
00076         }
00077     }
00078     else
00079     {
00080         forAll(points, pointI)
00081         {
00082             meshTools::writeOBJ(str, points[pointI]);
00083         }
00084 
00085         pointMap = identity(points.size());
00086     }
00087 
00088     forAll(edges, edgeI)
00089     {
00090         const edge& e = edges[edgeI];
00091 
00092         str<< "l " << pointMap[e[0]]+1 << ' ' << pointMap[e[1]]+1 << nl;
00093     }
00094 }
00095 
00096 
00097 //- Writes edges.
00098 void Foam::faceCoupleInfo::writeOBJ
00099 (
00100     const fileName& fName,
00101     const pointField& points0,
00102     const pointField& points1
00103 )
00104 {
00105     Pout<< "Writing connections as edges to " << fName << endl;
00106 
00107     OFstream str(fName);
00108 
00109     label vertI = 0;
00110 
00111     forAll(points0, i)
00112     {
00113         meshTools::writeOBJ(str, points0[i]);
00114         vertI++;
00115         meshTools::writeOBJ(str, points1[i]);
00116         vertI++;
00117         str << "l " << vertI-1 << ' ' << vertI << nl;
00118     }
00119 }
00120 
00121 
00122 //- Writes face and point connectivity as .obj files.
00123 void Foam::faceCoupleInfo::writePointsFaces() const
00124 {
00125     const indirectPrimitivePatch& m = masterPatch();
00126     const indirectPrimitivePatch& s = slavePatch();
00127     const primitiveFacePatch& c = cutFaces();
00128 
00129     // Patches
00130     {
00131         OFstream str("masterPatch.obj");
00132         Pout<< "Writing masterPatch to " << str.name() << endl;
00133         meshTools::writeOBJ(str, m.localFaces(), m.localPoints());
00134     }
00135     {
00136         OFstream str("slavePatch.obj");
00137         Pout<< "Writing slavePatch to " << str.name() << endl;
00138         meshTools::writeOBJ(str, s.localFaces(), s.localPoints());
00139     }
00140     {
00141         OFstream str("cutFaces.obj");
00142         Pout<< "Writing cutFaces to " << str.name() << endl;
00143         meshTools::writeOBJ(str, c.localFaces(), c.localPoints());
00144     }
00145 
00146     // Point connectivity
00147     {
00148         Pout<< "Writing cutToMasterPoints to cutToMasterPoints.obj" << endl;
00149 
00150         writeOBJ
00151         (
00152             "cutToMasterPoints.obj",
00153             m.localPoints(),
00154             pointField(c.localPoints(), masterToCutPoints_));
00155     }
00156     {
00157         Pout<< "Writing cutToSlavePoints to cutToSlavePoints.obj" << endl;
00158 
00159         writeOBJ
00160         (
00161             "cutToSlavePoints.obj",
00162             s.localPoints(),
00163             pointField(c.localPoints(), slaveToCutPoints_)
00164         );
00165     }
00166 
00167     // Face connectivity
00168     {
00169         Pout<< "Writing cutToMasterFaces to cutToMasterFaces.obj" << endl;
00170 
00171         pointField equivMasterFaces(c.size());
00172 
00173         forAll(cutToMasterFaces(), cutFaceI)
00174         {
00175             label masterFaceI = cutToMasterFaces()[cutFaceI];
00176 
00177             if (masterFaceI != -1)
00178             {
00179                 equivMasterFaces[cutFaceI] = m[masterFaceI].centre(m.points());
00180             }
00181             else
00182             {
00183                 WarningIn("writePointsFaces()")
00184                     << "No master face for cut face " << cutFaceI
00185                     << " at position " << c[cutFaceI].centre(c.points())
00186                     << endl;
00187 
00188                 equivMasterFaces[cutFaceI] = vector::zero;
00189             }
00190         }
00191 
00192         writeOBJ
00193         (
00194             "cutToMasterFaces.obj",
00195             calcFaceCentres<List>(c, cutPoints(), 0, c.size()),
00196             equivMasterFaces
00197         );
00198     }
00199 
00200     {
00201         Pout<< "Writing cutToSlaveFaces to cutToSlaveFaces.obj" << endl;
00202 
00203         pointField equivSlaveFaces(c.size());
00204 
00205         forAll(cutToSlaveFaces(), cutFaceI)
00206         {
00207             label slaveFaceI = cutToSlaveFaces()[cutFaceI];
00208 
00209             equivSlaveFaces[cutFaceI] = s[slaveFaceI].centre(s.points());
00210         }
00211 
00212         writeOBJ
00213         (
00214             "cutToSlaveFaces.obj",
00215             calcFaceCentres<List>(c, cutPoints(), 0, c.size()),
00216             equivSlaveFaces
00217         );
00218     }
00219 
00220     Pout<< endl;
00221 }
00222 
00223 
00224 void Foam::faceCoupleInfo::writeEdges
00225 (
00226     const labelList& cutToMasterEdges,
00227     const labelList& cutToSlaveEdges
00228 ) const
00229 {
00230     const indirectPrimitivePatch& m = masterPatch();
00231     const indirectPrimitivePatch& s = slavePatch();
00232     const primitiveFacePatch& c = cutFaces();
00233 
00234     // Edge connectivity
00235     {
00236         OFstream str("cutToMasterEdges.obj");
00237         Pout<< "Writing cutToMasterEdges to " << str.name() << endl;
00238 
00239         label vertI = 0;
00240 
00241         forAll(cutToMasterEdges, cutEdgeI)
00242         {
00243             if (cutToMasterEdges[cutEdgeI] != -1)
00244             {
00245                 const edge& masterEdge =
00246                     m.edges()[cutToMasterEdges[cutEdgeI]];
00247                 const edge& cutEdge = c.edges()[cutEdgeI];
00248 
00249                 meshTools::writeOBJ(str, m.localPoints()[masterEdge[0]]);
00250                 vertI++;
00251                 meshTools::writeOBJ(str, m.localPoints()[masterEdge[1]]);
00252                 vertI++;
00253                 meshTools::writeOBJ(str, c.localPoints()[cutEdge[0]]);
00254                 vertI++;
00255                 meshTools::writeOBJ(str, c.localPoints()[cutEdge[1]]);
00256                 vertI++;
00257                 str << "l " << vertI-3 << ' ' << vertI-2 << nl;
00258                 str << "l " << vertI-3 << ' ' << vertI-1 << nl;
00259                 str << "l " << vertI-3 << ' ' << vertI << nl;
00260                 str << "l " << vertI-2 << ' ' << vertI-1 << nl;
00261                 str << "l " << vertI-2 << ' ' << vertI << nl;
00262                 str << "l " << vertI-1 << ' ' << vertI << nl;
00263             }
00264         }
00265     }
00266     {
00267         OFstream str("cutToSlaveEdges.obj");
00268         Pout<< "Writing cutToSlaveEdges to " << str.name() << endl;
00269 
00270         label vertI = 0;
00271 
00272         labelList slaveToCut(invert(s.nEdges(), cutToSlaveEdges));
00273 
00274         forAll(slaveToCut, edgeI)
00275         {
00276             if (slaveToCut[edgeI] != -1)
00277             {
00278                 const edge& slaveEdge = s.edges()[edgeI];
00279                 const edge& cutEdge = c.edges()[slaveToCut[edgeI]];
00280 
00281                 meshTools::writeOBJ(str, s.localPoints()[slaveEdge[0]]);
00282                 vertI++;
00283                 meshTools::writeOBJ(str, s.localPoints()[slaveEdge[1]]);
00284                 vertI++;
00285                 meshTools::writeOBJ(str, c.localPoints()[cutEdge[0]]);
00286                 vertI++;
00287                 meshTools::writeOBJ(str, c.localPoints()[cutEdge[1]]);
00288                 vertI++;
00289                 str << "l " << vertI-3 << ' ' << vertI-2 << nl;
00290                 str << "l " << vertI-3 << ' ' << vertI-1 << nl;
00291                 str << "l " << vertI-3 << ' ' << vertI << nl;
00292                 str << "l " << vertI-2 << ' ' << vertI-1 << nl;
00293                 str << "l " << vertI-2 << ' ' << vertI << nl;
00294                 str << "l " << vertI-1 << ' ' << vertI << nl;
00295             }
00296         }
00297     }
00298 
00299     Pout<< endl;
00300 }
00301 
00302 
00303 // Given an edgelist and a map for the points on the edges it tries to find
00304 // the corresponding patch edges.
00305 Foam::labelList Foam::faceCoupleInfo::findMappedEdges
00306 (
00307     const edgeList& edges,
00308     const labelList& pointMap,
00309     const indirectPrimitivePatch& patch
00310 )
00311 {
00312     labelList toPatchEdges(edges.size());
00313 
00314     forAll(toPatchEdges, edgeI)
00315     {
00316         const edge& e = edges[edgeI];
00317 
00318         label v0 = pointMap[e[0]];
00319         label v1 = pointMap[e[1]];
00320 
00321         toPatchEdges[edgeI] =
00322             meshTools::findEdge
00323             (
00324                 patch.edges(),
00325                 patch.pointEdges()[v0],
00326                 v0,
00327                 v1
00328             );
00329     }
00330     return toPatchEdges;
00331 }
00332 
00333 
00334 // Detect a cut edge which originates from two boundary faces having different
00335 // polyPatches.
00336 bool Foam::faceCoupleInfo::regionEdge
00337 (
00338     const polyMesh& slaveMesh,
00339     const label slaveEdgeI
00340 ) const
00341 {
00342     const labelList& eFaces = slavePatch().edgeFaces()[slaveEdgeI];
00343 
00344     if (eFaces.size() == 1)
00345     {
00346         return true;
00347     }
00348     else
00349     {
00350         // Count how many different patches connected to this edge.
00351 
00352         label patch0 = -1;
00353 
00354         forAll(eFaces, i)
00355         {
00356             label faceI = eFaces[i];
00357 
00358             label meshFaceI = slavePatch().addressing()[faceI];
00359 
00360             label patchI = slaveMesh.boundaryMesh().whichPatch(meshFaceI);
00361 
00362             if (patch0 == -1)
00363             {
00364                 patch0 = patchI;
00365             }
00366             else if (patchI != patch0)
00367             {
00368                 // Found two different patches connected to this edge.
00369                 return true;
00370             }
00371         }
00372         return false;
00373     }
00374 }
00375 
00376 
00377 // Find edge using pointI that is most aligned with vector between
00378 // master points. Patchdivision tells us whether or not to use
00379 // patch information to match edges.
00380 Foam::label Foam::faceCoupleInfo::mostAlignedCutEdge
00381 (
00382     const bool report,
00383     const polyMesh& slaveMesh,
00384     const bool patchDivision,
00385     const labelList& cutToMasterEdges,
00386     const labelList& cutToSlaveEdges,
00387     const label pointI,
00388     const label edgeStart,
00389     const label edgeEnd
00390 ) const
00391 {
00392     const pointField& localPoints = cutFaces().localPoints();
00393 
00394     const labelList& pEdges = cutFaces().pointEdges()[pointI];
00395 
00396     if (report)
00397     {
00398         Pout<< "mostAlignedEdge : finding nearest edge among "
00399             << UIndirectList<edge>(cutFaces().edges(), pEdges)()
00400             << " connected to point " << pointI
00401             << " coord:" << localPoints[pointI]
00402             << " running between " << edgeStart << " coord:"
00403             << localPoints[edgeStart]
00404             << " and " << edgeEnd << " coord:"
00405             << localPoints[edgeEnd]
00406             << endl;
00407     }
00408 
00409     // Find the edge that gets us nearest end.
00410 
00411     label maxEdgeI = -1;
00412     scalar maxCos = -GREAT;
00413 
00414     forAll(pEdges, i)
00415     {
00416         label edgeI = pEdges[i];
00417 
00418         if
00419         (
00420            !(
00421                 patchDivision
00422              && cutToMasterEdges[edgeI] == -1
00423             )
00424          || (
00425                 patchDivision
00426              && regionEdge(slaveMesh, cutToSlaveEdges[edgeI])
00427             )
00428         )
00429         {
00430             const edge& e = cutFaces().edges()[edgeI];
00431 
00432             label otherPointI = e.otherVertex(pointI);
00433 
00434             if (otherPointI == edgeEnd)
00435             {
00436                 // Shortcut: found edge end point.
00437                 if (report)
00438                 {
00439                     Pout<< "    mostAlignedEdge : found end point " << edgeEnd
00440                         << endl;
00441                 }
00442                 return edgeI;
00443             }
00444 
00445             // Get angle between edge and edge to masterEnd
00446 
00447             vector eVec(localPoints[otherPointI] - localPoints[pointI]);
00448 
00449             scalar magEVec = mag(eVec);
00450 
00451             if (magEVec < VSMALL)
00452             {
00453                 WarningIn("faceCoupleInfo::mostAlignedEdge")
00454                     << "Crossing zero sized edge " << edgeI
00455                     << " coords:" << localPoints[otherPointI]
00456                     << localPoints[pointI]
00457                     << " when walking from " << localPoints[edgeStart]
00458                     << " to " << localPoints[edgeEnd]
00459                     << endl;
00460                 return edgeI;
00461             }
00462 
00463             eVec /= magEVec;
00464 
00465             vector eToEndPoint(localPoints[edgeEnd] - localPoints[otherPointI]);
00466             eToEndPoint /= mag(eToEndPoint);
00467 
00468             scalar cosAngle = eVec & eToEndPoint;
00469 
00470             if (report)
00471             {
00472                 Pout<< "    edge:" << e << " points:" << localPoints[pointI]
00473                     << localPoints[otherPointI]
00474                     << "  vec:" << eVec
00475                     << "  vecToEnd:" << eToEndPoint
00476                     << " cosAngle:" << cosAngle
00477                     << endl;
00478             }
00479 
00480             if (cosAngle > maxCos)
00481             {
00482                 maxCos = cosAngle;
00483                 maxEdgeI = edgeI;
00484             }
00485         }
00486     }
00487 
00488     if (maxCos > 1 - angleTol_)
00489     {
00490         return maxEdgeI;
00491     }
00492     else
00493     {
00494         return -1;
00495     }
00496 }
00497 
00498 
00499 // Construct points to split points map (in cut addressing)
00500 void Foam::faceCoupleInfo::setCutEdgeToPoints(const labelList& cutToMasterEdges)
00501 {
00502     labelListList masterToCutEdges
00503     (
00504         invertOneToMany
00505         (
00506             masterPatch().nEdges(),
00507             cutToMasterEdges
00508         )
00509     );
00510 
00511     const edgeList& cutEdges = cutFaces().edges();
00512 
00513     // Size extra big so searching is faster
00514     cutEdgeToPoints_.resize
00515     (
00516         masterPatch().nEdges()
00517       + slavePatch().nEdges()
00518       + cutEdges.size()
00519     );
00520 
00521     forAll(masterToCutEdges, masterEdgeI)
00522     {
00523         const edge& masterE = masterPatch().edges()[masterEdgeI];
00524 
00525         //Pout<< "Master:" << masterPatch().localPoints()[masterE[0]] << ' '
00526         //    << masterPatch().localPoints()[masterE[1]] << endl;
00527 
00528         const labelList& stringedEdges = masterToCutEdges[masterEdgeI];
00529 
00530         if (stringedEdges.empty())
00531         {
00532             FatalErrorIn
00533             (
00534                 "faceCoupleInfo::setCutEdgeToPoints"
00535                 "(const labelList&)"
00536             )   << "Did not match all of master edges to cutFace edges"
00537                 << nl
00538                 << "First unmatched edge:" << masterEdgeI << " endPoints:"
00539                 << masterPatch().localPoints()[masterE[0]]
00540                 << masterPatch().localPoints()[masterE[1]]
00541                 << endl
00542                 << "This usually means that the slave patch is not a"
00543                 << " subdivision of the master patch"
00544                 << abort(FatalError);
00545         }
00546         else if (stringedEdges.size() > 1)
00547         {
00548             // String up the edges between e[0] and e[1]. Store the points
00549             // inbetween e[0] and e[1] (all in cutFaces() labels)
00550 
00551             DynamicList<label> splitPoints(stringedEdges.size()-1);
00552 
00553             // Unsplit edge endpoints
00554             const edge unsplitEdge
00555             (
00556                 masterToCutPoints_[masterE[0]],
00557                 masterToCutPoints_[masterE[1]]
00558             );
00559 
00560             label startVertI = unsplitEdge[0];
00561             label startEdgeI = -1;
00562 
00563             while (startVertI != unsplitEdge[1])
00564             {
00565                 // Loop over all string of edges. Update
00566                 // - startVertI : previous vertex
00567                 // - startEdgeI : previous edge
00568                 // and insert any points into splitPoints
00569 
00570                 // For checking
00571                 label oldStart = startVertI;
00572 
00573                 forAll(stringedEdges, i)
00574                 {
00575                     label edgeI = stringedEdges[i];
00576 
00577                     if (edgeI != startEdgeI)
00578                     {
00579                         const edge& e = cutEdges[edgeI];
00580 
00581                         //Pout<< "    cut:" << e << " at:"
00582                         //    << cutFaces().localPoints()[e[0]]
00583                         //    << ' ' << cutFaces().localPoints()[e[1]] << endl;
00584 
00585                         if (e[0] == startVertI)
00586                         {
00587                             startEdgeI = edgeI;
00588                             startVertI = e[1];
00589                             if (e[1] != unsplitEdge[1])
00590                             {
00591                                 splitPoints.append(e[1]);
00592                             }
00593                             break;
00594                         }
00595                         else if (e[1] == startVertI)
00596                         {
00597                             startEdgeI = edgeI;
00598                             startVertI = e[0];
00599                             if (e[0] != unsplitEdge[1])
00600                             {
00601                                 splitPoints.append(e[0]);
00602                             }
00603                             break;
00604                         }
00605                     }
00606                 }
00607 
00608                 // Check
00609                 if (oldStart == startVertI)
00610                 {
00611                     FatalErrorIn
00612                     (
00613                         "faceCoupleInfo::setCutEdgeToPoints"
00614                         "(const labelList&)"
00615                     )   << " unsplitEdge:" << unsplitEdge
00616                         << " does not correspond to split edges "
00617                         << UIndirectList<edge>(cutEdges, stringedEdges)()
00618                         << abort(FatalError);
00619                 }
00620             }
00621 
00622             //Pout<< "For master edge:"
00623             //    << unsplitEdge
00624             //    << " Found stringed points "
00625             //    <<  UIndirectList<point>
00626             //        (
00627             //            cutFaces().localPoints(),
00628             //            splitPoints.shrink()
00629             //        )()
00630             //    << endl;
00631 
00632             cutEdgeToPoints_.insert(unsplitEdge, splitPoints.shrink());
00633         }
00634     }
00635 }
00636 
00637 
00638 // Determines rotation for f1 to match up with f0, i.e. the index in f0 of
00639 // the first point of f1.
00640 Foam::label Foam::faceCoupleInfo::matchFaces
00641 (
00642     const scalar absTol,
00643     const pointField& points0,
00644     const face& f0,
00645     const pointField& points1,
00646     const face& f1,
00647     const bool sameOrientation
00648 )
00649 {
00650     if (f0.size() != f1.size())
00651     {
00652         FatalErrorIn
00653         (
00654             "faceCoupleInfo::matchFaces"
00655             "(const scalar, const face&, const pointField&"
00656             ", const face&, const pointField&)"
00657         )   << "Different sizes for supposedly matching faces." << nl
00658             << "f0:" << f0 << " coords:" << UIndirectList<point>(points0, f0)()
00659             << nl
00660             << "f1:" << f1 << " coords:" << UIndirectList<point>(points1, f1)()
00661             << abort(FatalError);
00662     }
00663 
00664     const scalar absTolSqr = sqr(absTol);
00665 
00666 
00667     label matchFp = -1;
00668 
00669     forAll(f0, startFp)
00670     {
00671         // See -if starting from startFp on f0- the two faces match.
00672         bool fullMatch = true;
00673 
00674         label fp0 = startFp;
00675         label fp1 = 0;
00676 
00677         forAll(f1, i)
00678         {
00679             scalar distSqr = Foam::magSqr(points0[f0[fp0]] - points1[f1[fp1]]);
00680 
00681             if (distSqr > absTolSqr)
00682             {
00683                 fullMatch = false;
00684                 break;
00685             }
00686 
00687             fp0 = f0.fcIndex(fp0);  // walk forward
00688 
00689             if (sameOrientation)
00690             {
00691                 fp1 = f1.fcIndex(fp1);
00692             }
00693             else
00694             {
00695                 fp1 = f1.rcIndex(fp1);
00696             }
00697         }
00698 
00699         if (fullMatch)
00700         {
00701             matchFp = startFp;
00702             break;
00703         }
00704     }
00705 
00706     if (matchFp == -1)
00707     {
00708         FatalErrorIn
00709         (
00710             "faceCoupleInfo::matchFaces"
00711             "(const scalar, const face&, const pointField&"
00712             ", const face&, const pointField&)"
00713         )   << "No unique match between two faces" << nl
00714             << "Face " << f0 << " coords "
00715             << UIndirectList<point>(points0, f0)() << nl
00716             << "Face " << f1 << " coords "
00717             << UIndirectList<point>(points1, f1)()
00718             << "when using tolerance " << absTol
00719             << " and forwardMatching:" << sameOrientation
00720             << abort(FatalError);
00721     }
00722 
00723     return matchFp;
00724 }
00725 
00726 
00727 // Find correspondence from patch points to cut points. This might
00728 // detect shared points so the output is a patch-to-cut point list
00729 // and a compaction list for the cut points (which will always be equal or more
00730 // connected than the patch).
00731 // Returns true if there are any duplicates.
00732 bool Foam::faceCoupleInfo::matchPointsThroughFaces
00733 (
00734     const scalar absTol,
00735     const pointField& cutPoints,
00736     const faceList& cutFaces,
00737     const pointField& patchPoints,
00738     const faceList& patchFaces,
00739     const bool sameOrientation,
00740 
00741     labelList& patchToCutPoints,    // patch to (uncompacted) cut points
00742     labelList& cutToCompact,        // compaction list for cut points
00743     labelList& compactToCut         // inverse ,,
00744 )
00745 {
00746 
00747     // From slave to cut point
00748     patchToCutPoints.setSize(patchPoints.size());
00749     patchToCutPoints = -1;
00750 
00751     // Compaction list for cut points: either -1 or index into master which
00752     // gives the point to compact to.
00753     labelList cutPointRegion(cutPoints.size(), -1);
00754     DynamicList<label> cutPointRegionMaster;
00755 
00756     forAll(patchFaces, patchFaceI)
00757     {
00758         const face& patchF = patchFaces[patchFaceI];
00759 
00760         //const face& cutF = cutFaces[patchToCutFaces[patchFaceI]];
00761         const face& cutF = cutFaces[patchFaceI];
00762 
00763         // Do geometric matching to get position of cutF[0] in patchF
00764         label patchFp = matchFaces
00765         (
00766             absTol,
00767             patchPoints,
00768             patchF,
00769             cutPoints,
00770             cutF,
00771             sameOrientation        // orientation
00772         );
00773 
00774         forAll(cutF, cutFp)
00775         {
00776             label cutPointI = cutF[cutFp];
00777             label patchPointI = patchF[patchFp];
00778 
00779             //const point& cutPt = cutPoints[cutPointI];
00780             //const point& patchPt = patchPoints[patchPointI];
00781             //if (mag(cutPt - patchPt) > SMALL)
00782             //{
00783             //    FatalErrorIn("matchPointsThroughFaces")
00784             //    << "cutP:" << cutPt
00785             //    << " patchP:" << patchPt
00786             //    << abort(FatalError);
00787             //}
00788 
00789             if (patchToCutPoints[patchPointI] == -1)
00790             {
00791                 patchToCutPoints[patchPointI] = cutPointI;
00792             }
00793             else if (patchToCutPoints[patchPointI] != cutPointI)
00794             {
00795                 // Multiple cut points connecting to same patch.
00796                 // Check if already have region & region master for this set
00797                 label otherCutPointI = patchToCutPoints[patchPointI];
00798 
00799                 //Pout<< "PatchPoint:" << patchPt
00800                 //    << " matches to:" << cutPointI
00801                 //    << " coord:" << cutPoints[cutPointI]
00802                 //    << " and to:" << otherCutPointI
00803                 //    << " coord:" << cutPoints[otherCutPointI]
00804                 //    << endl;
00805 
00806                 if (cutPointRegion[otherCutPointI] != -1)
00807                 {
00808                     // Have region for this set. Copy.
00809                     label region = cutPointRegion[otherCutPointI];
00810                     cutPointRegion[cutPointI] = region;
00811 
00812                     // Update region master with min point label
00813                     cutPointRegionMaster[region] = min
00814                     (
00815                         cutPointRegionMaster[region],
00816                         cutPointI
00817                     );
00818                 }
00819                 else
00820                 {
00821                     // Create new region.
00822                     label region = cutPointRegionMaster.size();
00823                     cutPointRegionMaster.append
00824                     (
00825                         min(cutPointI, otherCutPointI)
00826                     );
00827                     cutPointRegion[cutPointI] = region;
00828                     cutPointRegion[otherCutPointI] = region;
00829                 }
00830             }
00831 
00832             if (sameOrientation)
00833             {
00834                 patchFp = patchF.fcIndex(patchFp);
00835             }
00836             else
00837             {
00838                 patchFp = patchF.rcIndex(patchFp);
00839             }
00840         }
00841     }
00842 
00843     // Rework region&master into compaction array
00844     compactToCut.setSize(cutPointRegion.size());
00845     cutToCompact.setSize(cutPointRegion.size());
00846     cutToCompact = -1;
00847     label compactPointI = 0;
00848 
00849     forAll(cutPointRegion, i)
00850     {
00851         if (cutPointRegion[i] == -1)
00852         {
00853             // Unduplicated point. Allocate new compacted point.
00854             cutToCompact[i] = compactPointI;
00855             compactToCut[compactPointI] = i;
00856             compactPointI++;
00857         }
00858         else
00859         {
00860             // Duplicate point. Get master.
00861 
00862             label masterPointI = cutPointRegionMaster[cutPointRegion[i]];
00863 
00864             if (cutToCompact[masterPointI] == -1)
00865             {
00866                 cutToCompact[masterPointI] = compactPointI;
00867                 compactToCut[compactPointI] = masterPointI;
00868                 compactPointI++;
00869             }
00870             cutToCompact[i] = cutToCompact[masterPointI];
00871         }
00872     }
00873     compactToCut.setSize(compactPointI);
00874 
00875     return compactToCut.size() != cutToCompact.size();
00876 }
00877 
00878 
00879 // Return max distance from any point on cutF to masterF
00880 Foam::scalar Foam::faceCoupleInfo::maxDistance
00881 (
00882     const face& cutF,
00883     const pointField& cutPoints,
00884     const face& masterF,
00885     const pointField& masterPoints
00886 )
00887 {
00888     scalar maxDist = -GREAT;
00889 
00890     forAll(cutF, fp)
00891     {
00892         const point& cutPt = cutPoints[cutF[fp]];
00893 
00894         pointHit pHit = masterF.nearestPoint(cutPt, masterPoints);
00895 
00896         maxDist = max(maxDist, pHit.distance());
00897     }
00898     return maxDist;
00899 }
00900 
00901 
00902 void Foam::faceCoupleInfo::findPerfectMatchingFaces
00903 (
00904     const primitiveMesh& mesh0,
00905     const primitiveMesh& mesh1,
00906     const scalar absTol,
00907 
00908     labelList& mesh0Faces,
00909     labelList& mesh1Faces
00910 )
00911 {
00912     // Face centres of external faces (without invoking
00913     // mesh.faceCentres since mesh might have been clearedOut)
00914 
00915     pointField fc0
00916     (
00917         calcFaceCentres<List>
00918         (
00919             mesh0.faces(),
00920             mesh0.points(),
00921             mesh0.nInternalFaces(),
00922             mesh0.nFaces() - mesh0.nInternalFaces()
00923         )
00924     );
00925 
00926     pointField fc1
00927     (
00928         calcFaceCentres<List>
00929         (
00930             mesh1.faces(),
00931             mesh1.points(),
00932             mesh1.nInternalFaces(),
00933             mesh1.nFaces() - mesh1.nInternalFaces()
00934         )
00935     );
00936 
00937 
00938     if (debug)
00939     {
00940         Pout<< "Face matching tolerance : " << absTol << endl;
00941     }
00942 
00943 
00944     // Match geometrically
00945     labelList from1To0;
00946     bool matchedAllFaces = matchPoints
00947     (
00948         fc1,
00949         fc0,
00950         scalarField(fc1.size(), absTol),
00951         false,
00952         from1To0
00953     );
00954 
00955     if (matchedAllFaces)
00956     {
00957         Warning
00958             << "faceCoupleInfo::faceCoupleInfo : "
00959             << "Matched ALL " << fc1.size()
00960             << " boundary faces of mesh0 to boundary faces of mesh1." << endl
00961             << "This is only valid if the mesh to add is fully"
00962             << " enclosed by the mesh it is added to." << endl;
00963     }
00964 
00965 
00966     // Collect matches.
00967     label nMatched = 0;
00968 
00969     mesh0Faces.setSize(fc0.size());
00970     mesh1Faces.setSize(fc1.size());
00971 
00972     forAll(from1To0, i)
00973     {
00974         if (from1To0[i] != -1)
00975         {
00976             mesh1Faces[nMatched] = i + mesh1.nInternalFaces();
00977             mesh0Faces[nMatched] = from1To0[i] + mesh0.nInternalFaces();
00978 
00979             nMatched++;
00980         }
00981     }
00982 
00983     mesh0Faces.setSize(nMatched);
00984     mesh1Faces.setSize(nMatched);
00985 }
00986 
00987 
00988 void Foam::faceCoupleInfo::findSlavesCoveringMaster
00989 (
00990     const primitiveMesh& mesh0,
00991     const primitiveMesh& mesh1,
00992     const scalar absTol,
00993 
00994     labelList& mesh0Faces,
00995     labelList& mesh1Faces
00996 )
00997 {
00998     // Construct octree from all mesh0 boundary faces
00999     labelList bndFaces(mesh0.nFaces()-mesh0.nInternalFaces());
01000     forAll(bndFaces, i)
01001     {
01002         bndFaces[i] = mesh0.nInternalFaces() + i;
01003     }
01004 
01005     treeBoundBox overallBb(mesh0.points());
01006 
01007     Random rndGen(123456);
01008 
01009     indexedOctree<treeDataFace> tree
01010     (
01011         treeDataFace    // all information needed to search faces
01012         (
01013             false,                      // do not cache bb
01014             mesh0,
01015             bndFaces                    // boundary faces only
01016         ),
01017         overallBb.extend(rndGen, 1E-4), // overall search domain
01018         8,                              // maxLevel
01019         10,                             // leafsize
01020         3.0                             // duplicity
01021     );
01022 
01023     if (debug)
01024     {
01025         Pout<< "findSlavesCoveringMaster :"
01026             << " constructed octree for mesh0 boundary faces" << endl;
01027     }
01028 
01029     // Search nearest face for every mesh1 boundary face.
01030 
01031     labelHashSet mesh0Set(mesh0.nFaces() - mesh0.nInternalFaces());
01032     labelHashSet mesh1Set(mesh1.nFaces() - mesh1.nInternalFaces());
01033 
01034     for
01035     (
01036         label mesh1FaceI = mesh1.nInternalFaces();
01037         mesh1FaceI < mesh1.nFaces();
01038         mesh1FaceI++
01039     )
01040     {
01041         const face& f1 = mesh1.faces()[mesh1FaceI];
01042 
01043         // Generate face centre (prevent cellCentres() reconstruction)
01044         point fc(f1.centre(mesh1.points()));
01045 
01046         pointIndexHit nearInfo = tree.findNearest(fc, Foam::sqr(absTol));
01047 
01048         if (nearInfo.hit())
01049         {
01050             label mesh0FaceI = tree.shapes().faceLabels()[nearInfo.index()];
01051 
01052             // Check if points of f1 actually lie on top of mesh0 face
01053             // This is the bit that might fail since if f0 is severely warped
01054             // and f1's points are not present on f0 (since subdivision)
01055             // then the distance of the points to f0 might be quite large
01056             // and the test will fail. This all should in fact be some kind
01057             // of walk from known corresponding points/edges/faces.
01058             scalar dist =
01059                 maxDistance
01060                 (
01061                     f1,
01062                     mesh1.points(),
01063                     mesh0.faces()[mesh0FaceI],
01064                     mesh0.points()
01065                 );
01066 
01067             if (dist < absTol)
01068             {
01069                 mesh0Set.insert(mesh0FaceI);
01070                 mesh1Set.insert(mesh1FaceI);
01071             }
01072         }
01073     }
01074 
01075     if (debug)
01076     {
01077         Pout<< "findSlavesCoveringMaster :"
01078             << " matched " << mesh1Set.size() << " mesh1 faces to "
01079             << mesh0Set.size() << " mesh0 faces" << endl;
01080     }
01081 
01082     mesh0Faces = mesh0Set.toc();
01083     mesh1Faces = mesh1Set.toc();
01084 }
01085 
01086 
01087 // Grow cutToMasterFace across 'internal' edges.
01088 Foam::label Foam::faceCoupleInfo::growCutFaces
01089 (
01090     const labelList& cutToMasterEdges,
01091     Map<labelList>& candidates
01092 )
01093 {
01094     if (debug)
01095     {
01096         Pout<< "growCutFaces :"
01097             << " growing cut faces to masterPatch" << endl;
01098     }
01099 
01100     label nTotChanged = 0;
01101 
01102     while (true)
01103     {
01104         const labelListList& cutFaceEdges = cutFaces().faceEdges();
01105         const labelListList& cutEdgeFaces = cutFaces().edgeFaces();
01106 
01107         label nChanged = 0;
01108 
01109         forAll(cutToMasterFaces_, cutFaceI)
01110         {
01111             const label masterFaceI = cutToMasterFaces_[cutFaceI];
01112 
01113             if (masterFaceI != -1)
01114             {
01115                 // We now have a cutFace for which we have already found a
01116                 // master face. Grow this masterface across any internal edge
01117                 // (internal: no corresponding master edge)
01118 
01119                 const labelList& fEdges = cutFaceEdges[cutFaceI];
01120 
01121                 forAll(fEdges, i)
01122                 {
01123                     const label cutEdgeI = fEdges[i];
01124 
01125                     if (cutToMasterEdges[cutEdgeI] == -1)
01126                     {
01127                         // So this edge:
01128                         // - internal to the cutPatch (since all region edges
01129                         //   marked before)
01130                         // - on cutFaceI which corresponds to masterFace.
01131                         // Mark all connected faces with this masterFace.
01132 
01133                         const labelList& eFaces = cutEdgeFaces[cutEdgeI];
01134 
01135                         forAll(eFaces, j)
01136                         {
01137                             const label faceI = eFaces[j];
01138 
01139                             if (cutToMasterFaces_[faceI] == -1)
01140                             {
01141                                 cutToMasterFaces_[faceI] = masterFaceI;
01142                                 candidates.erase(faceI);
01143                                 nChanged++;
01144                             }
01145                             else if (cutToMasterFaces_[faceI] != masterFaceI)
01146                             {
01147                                 const pointField& cutPoints =
01148                                     cutFaces().points();
01149                                 const pointField& masterPoints =
01150                                     masterPatch().points();
01151 
01152                                 const edge& e = cutFaces().edges()[cutEdgeI];
01153 
01154                                 label myMaster = cutToMasterFaces_[faceI];
01155                                 const face& myF = masterPatch()[myMaster];
01156 
01157                                 const face& nbrF = masterPatch()[masterFaceI];
01158 
01159                                 FatalErrorIn
01160                                 (
01161                                     "faceCoupleInfo::growCutFaces"
01162                                     "(const labelList&, Map<labelList>&)"
01163                                 )   << "Cut face "
01164                                     << cutFaces()[faceI].points(cutPoints)
01165                                     << " has master "
01166                                     << myMaster
01167                                     << " but also connects to nbr face "
01168                                     << cutFaces()[cutFaceI].points(cutPoints)
01169                                     << " with master " << masterFaceI
01170                                     << nl
01171                                     << "myMasterFace:"
01172                                     << myF.points(masterPoints)
01173                                     << "  nbrMasterFace:"
01174                                     << nbrF.points(masterPoints) << nl
01175                                     << "Across cut edge " << cutEdgeI
01176                                     << " coord:"
01177                                     << cutFaces().localPoints()[e[0]]
01178                                     << cutFaces().localPoints()[e[1]]
01179                                     << abort(FatalError);
01180                             }
01181                         }
01182                     }
01183                 }
01184             }
01185         }
01186 
01187         if (debug)
01188         {
01189             Pout<< "growCutFaces : Grown an additional " << nChanged
01190                 << " cut to master face" << " correspondence" << endl;
01191         }
01192 
01193         nTotChanged += nChanged;
01194 
01195         if (nChanged == 0)
01196         {
01197             break;
01198         }
01199     }
01200 
01201     return nTotChanged;
01202 }
01203 
01204 
01205 void Foam::faceCoupleInfo::checkMatch(const labelList& cutToMasterEdges) const
01206 {
01207     const pointField& cutLocalPoints = cutFaces().localPoints();
01208 
01209     const pointField& masterLocalPoints = masterPatch().localPoints();
01210     const faceList& masterLocalFaces = masterPatch().localFaces();
01211 
01212     forAll(cutToMasterEdges, cutEdgeI)
01213     {
01214         const edge& e = cutFaces().edges()[cutEdgeI];
01215 
01216         if (cutToMasterEdges[cutEdgeI] == -1)
01217         {
01218             // Internal edge. Check that master face is same on both sides.
01219             const labelList& cutEFaces = cutFaces().edgeFaces()[cutEdgeI];
01220 
01221             label masterFaceI = -1;
01222 
01223             forAll(cutEFaces, i)
01224             {
01225                 label cutFaceI = cutEFaces[i];
01226 
01227                 if (cutToMasterFaces_[cutFaceI] != -1)
01228                 {
01229                     if (masterFaceI == -1)
01230                     {
01231                         masterFaceI = cutToMasterFaces_[cutFaceI];
01232                     }
01233                     else if (masterFaceI != cutToMasterFaces_[cutFaceI])
01234                     {
01235                         label myMaster = cutToMasterFaces_[cutFaceI];
01236                         const face& myF = masterLocalFaces[myMaster];
01237 
01238                         const face& nbrF = masterLocalFaces[masterFaceI];
01239 
01240                         FatalErrorIn
01241                         (
01242                             "faceCoupleInfo::checkMatch(const labelList&) const"
01243                         )
01244                             << "Internal CutEdge " << e
01245                             << " coord:"
01246                             << cutLocalPoints[e[0]]
01247                             << cutLocalPoints[e[1]]
01248                             << " connects to master " << myMaster
01249                             << " and to master " << masterFaceI << nl
01250                             << "myMasterFace:"
01251                             << myF.points(masterLocalPoints)
01252                             << "  nbrMasterFace:"
01253                             << nbrF.points(masterLocalPoints)
01254                             << abort(FatalError);
01255                     }
01256                 }
01257             }
01258         }
01259     }
01260 }
01261 
01262 
01263 // Extends matching information by elimination across cutFaces using more
01264 // than one region edge. Updates cutToMasterFaces_ and sets candidates
01265 // which is for every cutface on a region edge the possible master faces.
01266 Foam::label Foam::faceCoupleInfo::matchEdgeFaces
01267 (
01268     const labelList& cutToMasterEdges,
01269     Map<labelList>& candidates
01270 )
01271 {
01272     // For every unassigned cutFaceI the possible list of master faces.
01273     candidates.clear();
01274     candidates.resize(cutFaces().size());
01275 
01276     label nChanged = 0;
01277 
01278     forAll(cutToMasterEdges, cutEdgeI)
01279     {
01280         label masterEdgeI = cutToMasterEdges[cutEdgeI];
01281 
01282         if (masterEdgeI != -1)
01283         {
01284             // So cutEdgeI is matched to masterEdgeI. For all cut faces using
01285             // the cut edge store the master face as a candidate match.
01286 
01287             const labelList& cutEFaces = cutFaces().edgeFaces()[cutEdgeI];
01288             const labelList& masterEFaces =
01289                 masterPatch().edgeFaces()[masterEdgeI];
01290 
01291             forAll(cutEFaces, i)
01292             {
01293                 label cutFaceI = cutEFaces[i];
01294 
01295                 if (cutToMasterFaces_[cutFaceI] == -1)
01296                 {
01297                     // So this face (cutFaceI) is on an edge (cutEdgeI) that
01298                     // is used by some master faces (masterEFaces). Check if
01299                     // the cutFaceI and some of these masterEFaces share more
01300                     // than one edge (which uniquely defines face)
01301 
01302                     // Combine master faces with current set of candidate
01303                     // master faces.
01304                     Map<labelList>::iterator fnd = candidates.find(cutFaceI);
01305 
01306                     if (fnd == candidates.end())
01307                     {
01308                         // No info yet for cutFaceI. Add all master faces as
01309                         // candidates
01310                         candidates.insert(cutFaceI, masterEFaces);
01311                     }
01312                     else
01313                     {
01314                         // From some other cutEdgeI there are already some
01315                         // candidate master faces. Check the overlap with
01316                         // the current set of master faces.
01317                         const labelList& masterFaces = fnd();
01318 
01319                         DynamicList<label> newCandidates(masterFaces.size());
01320 
01321                         forAll(masterEFaces, j)
01322                         {
01323                             if (findIndex(masterFaces, masterEFaces[j]) != -1)
01324                             {
01325                                 newCandidates.append(masterEFaces[j]);
01326                             }
01327                         }
01328 
01329                         if (newCandidates.size() == 1)
01330                         {
01331                             // We found a perfect match. Delete entry from
01332                             // candidates map.
01333                             cutToMasterFaces_[cutFaceI] = newCandidates[0];
01334                             candidates.erase(cutFaceI);
01335                             nChanged++;
01336                         }
01337                         else
01338                         {
01339                             // Should not happen?
01340                             //Pout<< "On edge:" << cutEdgeI
01341                             //    << " have connected masterFaces:"
01342                             //    << masterEFaces
01343                             //    << " and from previous edge we have"
01344                             //    << " connected masterFaces" << masterFaces
01345                             //    << " . Overlap:" << newCandidates << endl;
01346 
01347                             fnd() = newCandidates.shrink();
01348                         }
01349                     }
01350                 }
01351 
01352             }
01353         }
01354     }
01355 
01356     if (debug)
01357     {
01358         Pout<< "matchEdgeFaces : Found " << nChanged
01359             << " faces where there was"
01360             << " only one remaining choice for cut-master correspondence"
01361             << endl;
01362     }
01363 
01364     return nChanged;
01365 }
01366 
01367 
01368 // Gets a list of cutFaces (that use a master edge) and the candidate
01369 // master faces.
01370 // Finds most aligned master face.
01371 Foam::label Foam::faceCoupleInfo::geometricMatchEdgeFaces
01372 (
01373     Map<labelList>& candidates
01374 )
01375 {
01376     const pointField& cutPoints = cutFaces().points();
01377 
01378     label nChanged = 0;
01379 
01380     // Mark all master faces that have been matched so far.
01381 
01382     labelListList masterToCutFaces
01383     (
01384         invertOneToMany
01385         (
01386             masterPatch().size(),
01387             cutToMasterFaces_
01388         )
01389     );
01390 
01391     forAllConstIter(Map<labelList>, candidates, iter)
01392     {
01393         label cutFaceI = iter.key();
01394 
01395         const face& cutF = cutFaces()[cutFaceI];
01396 
01397         if (cutToMasterFaces_[cutFaceI] == -1)
01398         {
01399             const labelList& masterFaces = iter();
01400 
01401             // Find the best matching master face.
01402             scalar minDist = GREAT;
01403             label minMasterFaceI = -1;
01404 
01405             forAll(masterFaces, i)
01406             {
01407                 label masterFaceI = masterFaces[i];
01408 
01409                 if (masterToCutFaces[masterFaces[i]].empty())
01410                 {
01411                     scalar dist = maxDistance
01412                     (
01413                         cutF,
01414                         cutPoints,
01415                         masterPatch()[masterFaceI],
01416                         masterPatch().points()
01417                     );
01418 
01419                     if (dist < minDist)
01420                     {
01421                         minDist = dist;
01422                         minMasterFaceI = masterFaceI;
01423                     }
01424                 }
01425             }
01426 
01427             if (minMasterFaceI != -1)
01428             {
01429                 cutToMasterFaces_[cutFaceI] = minMasterFaceI;
01430                 masterToCutFaces[minMasterFaceI] = cutFaceI;
01431                 nChanged++;
01432             }
01433         }
01434     }
01435 
01436     // (inefficiently) force candidates to be uptodate.
01437     forAll(cutToMasterFaces_, cutFaceI)
01438     {
01439         if (cutToMasterFaces_[cutFaceI] != -1)
01440         {
01441             candidates.erase(cutFaceI);
01442         }
01443     }
01444 
01445 
01446     if (debug)
01447     {
01448         Pout<< "geometricMatchEdgeFaces : Found " << nChanged
01449             << " faces where there was"
01450             << " only one remaining choice for cut-master correspondence"
01451             << endl;
01452     }
01453 
01454     return nChanged;
01455 }
01456 
01457 
01458 // Calculate the set of cut faces inbetween master and slave patch
01459 // assuming perfect match (and optional face ordering on slave)
01460 void Foam::faceCoupleInfo::perfectPointMatch
01461 (
01462     const scalar absTol,
01463     const bool slaveFacesOrdered
01464 )
01465 {
01466     if (debug)
01467     {
01468         Pout<< "perfectPointMatch :"
01469             << " Matching master and slave to cut."
01470             << " Master and slave faces are identical" << nl;
01471 
01472         if (slaveFacesOrdered)
01473         {
01474             Pout<< "and master and slave faces are ordered"
01475                 << " (on coupled patches)" << endl;
01476         }
01477         else
01478         {
01479             Pout<< "and master and slave faces are not ordered"
01480                 << " (on coupled patches)" << endl;
01481         }
01482     }
01483 
01484     cutToMasterFaces_ = identity(masterPatch().size());
01485     cutPoints_ = masterPatch().localPoints();
01486     cutFacesPtr_.reset
01487     (
01488         new primitiveFacePatch
01489         (
01490             masterPatch().localFaces(),
01491             cutPoints_
01492         )
01493     );
01494     masterToCutPoints_ = identity(cutPoints_.size());
01495 
01496 
01497     // Cut faces to slave patch.
01498     bool matchedAllFaces = false;
01499 
01500     if (slaveFacesOrdered)
01501     {
01502         cutToSlaveFaces_ = identity(cutFaces().size());
01503         matchedAllFaces = (cutFaces().size() == slavePatch().size());
01504     }
01505     else
01506     {
01507         // Faces do not have to be ordered (but all have
01508         // to match). Note: Faces will be already ordered if we enter here from
01509         // construct from meshes.
01510         matchedAllFaces = matchPoints
01511         (
01512             calcFaceCentres<List>
01513             (
01514                 cutFaces(),
01515                 cutPoints_,
01516                 0,
01517                 cutFaces().size()
01518             ),
01519             calcFaceCentres<IndirectList>
01520             (
01521                 slavePatch(),
01522                 slavePatch().points(),
01523                 0,
01524                 slavePatch().size()
01525             ),
01526             scalarField(slavePatch().size(), absTol),
01527             true,
01528             cutToSlaveFaces_
01529         );
01530     }
01531 
01532 
01533     if (!matchedAllFaces)
01534     {
01535         FatalErrorIn
01536         (
01537             "faceCoupleInfo::perfectPointMatch"
01538             "(const scalar&, const bool)"
01539         )   << "Did not match all of the master faces to the slave faces"
01540             << endl
01541             << "This usually means that the slave patch and master patch"
01542             << " do not align to within " << absTol << " meter."
01543             << abort(FatalError);
01544     }
01545 
01546 
01547     // Find correspondence from slave points to cut points. This might
01548     // detect shared points so the output is a slave-to-cut point list
01549     // and a compaction list.
01550 
01551     labelList cutToCompact, compactToCut;
01552     matchPointsThroughFaces
01553     (
01554         absTol,
01555         cutFaces().localPoints(),
01556         reorder(cutToSlaveFaces_, cutFaces().localFaces()),
01557         slavePatch().localPoints(),
01558         slavePatch().localFaces(),
01559         false,                      // slave and cut have opposite orientation
01560 
01561         slaveToCutPoints_,          // slave to (uncompacted) cut points
01562         cutToCompact,               // compaction map: from cut to compacted
01563         compactToCut                // compaction map: from compacted to cut
01564     );
01565 
01566 
01567     // Use compaction lists to renumber cutPoints.
01568     cutPoints_ = UIndirectList<point>(cutPoints_, compactToCut)();
01569     {
01570         const faceList& cutLocalFaces = cutFaces().localFaces();
01571 
01572         faceList compactFaces(cutLocalFaces.size());
01573         forAll(cutLocalFaces, i)
01574         {
01575             compactFaces[i] = renumber(cutToCompact, cutLocalFaces[i]);
01576         }
01577         cutFacesPtr_.reset
01578         (
01579             new primitiveFacePatch
01580             (
01581                 compactFaces,
01582                 cutPoints_
01583             )
01584         );
01585     }
01586     inplaceRenumber(cutToCompact, slaveToCutPoints_);
01587     inplaceRenumber(cutToCompact, masterToCutPoints_);
01588 }
01589 
01590 
01591 // Calculate the set of cut faces inbetween master and slave patch
01592 // assuming that slave patch is subdivision of masterPatch.
01593 void Foam::faceCoupleInfo::subDivisionMatch
01594 (
01595     const polyMesh& slaveMesh,
01596     const bool patchDivision,
01597     const scalar absTol
01598 )
01599 {
01600     if (debug)
01601     {
01602         Pout<< "subDivisionMatch :"
01603             << " Matching master and slave to cut."
01604             << " Slave can be subdivision of master but all master edges"
01605             << " have to be covered by slave edges." << endl;
01606     }
01607 
01608     // Assume slave patch is subdivision of the master patch so cutFaces is a
01609     // copy of the slave (but reversed since orientation has to be like master).
01610     cutPoints_ = slavePatch().localPoints();
01611     {
01612         faceList cutFaces(slavePatch().size());
01613 
01614         forAll(cutFaces, i)
01615         {
01616             cutFaces[i] = slavePatch().localFaces()[i].reverseFace();
01617         }
01618         cutFacesPtr_.reset(new primitiveFacePatch(cutFaces, cutPoints_));
01619     }
01620 
01621     // Cut is copy of slave so addressing to slave is trivial.
01622     cutToSlaveFaces_ = identity(cutFaces().size());
01623     slaveToCutPoints_ = identity(slavePatch().nPoints());
01624 
01625     // Cut edges to slave patch
01626     labelList cutToSlaveEdges
01627     (
01628         findMappedEdges
01629         (
01630             cutFaces().edges(),
01631             slaveToCutPoints_,  //note:should be cutToSlavePoints but since iden
01632             slavePatch()
01633         )
01634     );
01635 
01636 
01637     if (debug)
01638     {
01639         OFstream str("regionEdges.obj");
01640 
01641         Pout<< "subDivisionMatch :"
01642             << " addressing for slave patch fully done."
01643             << " Dumping region edges to " << str.name() << endl;
01644 
01645         label vertI = 0;
01646 
01647         forAll(slavePatch().edges(), slaveEdgeI)
01648         {
01649             if (regionEdge(slaveMesh, slaveEdgeI))
01650             {
01651                 const edge& e = slavePatch().edges()[slaveEdgeI];
01652 
01653                 meshTools::writeOBJ(str, slavePatch().localPoints()[e[0]]);
01654                 vertI++;
01655                 meshTools::writeOBJ(str, slavePatch().localPoints()[e[1]]);
01656                 vertI++;
01657                 str<< "l " << vertI-1 << ' ' << vertI << nl;
01658             }
01659         }
01660     }
01661 
01662 
01663     // Addressing from cut to master.
01664 
01665     // Cut points to master patch
01666     // - determine master points to cut points. Has to be full!
01667     // - invert to get cut to master
01668     if (debug)
01669     {
01670         Pout<< "subDivisionMatch :"
01671             << " matching master points to cut points" << endl;
01672     }
01673 
01674     bool matchedAllPoints = matchPoints
01675     (
01676         masterPatch().localPoints(),
01677         cutPoints_,
01678         scalarField(masterPatch().nPoints(), absTol),
01679         true,
01680         masterToCutPoints_
01681     );
01682 
01683     if (!matchedAllPoints)
01684     {
01685         FatalErrorIn
01686         (
01687             "faceCoupleInfo::subDivisionMatch"
01688             "(const polyMesh&, const bool, const scalar)"
01689         )   << "Did not match all of the master points to the slave points"
01690             << endl
01691             << "This usually means that the slave patch is not a"
01692             << " subdivision of the master patch"
01693             << abort(FatalError);
01694     }
01695 
01696 
01697     // Do masterEdges to cutEdges :
01698     // - mark all edges between two masterEdge endpoints. (geometric test since
01699     //   this is the only distinction)
01700     // - this gives the borders inbetween which all cutfaces come from
01701     //   a single master face.
01702     if (debug)
01703     {
01704         Pout<< "subDivisionMatch :"
01705             << " matching cut edges to master edges" << endl;
01706     }
01707 
01708     const edgeList& masterEdges = masterPatch().edges();
01709     const pointField& masterPoints = masterPatch().localPoints();
01710 
01711     const edgeList& cutEdges = cutFaces().edges();
01712 
01713     labelList cutToMasterEdges(cutFaces().nEdges(), -1);
01714 
01715     forAll(masterEdges, masterEdgeI)
01716     {
01717         const edge& masterEdge = masterEdges[masterEdgeI];
01718 
01719         label cutPoint0 = masterToCutPoints_[masterEdge[0]];
01720         label cutPoint1 = masterToCutPoints_[masterEdge[1]];
01721 
01722         // Find edges between cutPoint0 and cutPoint1.
01723 
01724         label cutPointI = cutPoint0;
01725         do
01726         {
01727             // Find edge (starting at pointI on cut), aligned with master
01728             // edge.
01729             label cutEdgeI =
01730                 mostAlignedCutEdge
01731                 (
01732                     false,
01733                     slaveMesh,
01734                     patchDivision,
01735                     cutToMasterEdges,
01736                     cutToSlaveEdges,
01737                     cutPointI,
01738                     cutPoint0,
01739                     cutPoint1
01740                 );
01741 
01742             if (cutEdgeI == -1)
01743             {
01744                 // Problem. Report while matching to produce nice error message
01745                 mostAlignedCutEdge
01746                 (
01747                     true,
01748                     slaveMesh,
01749                     patchDivision,
01750                     cutToMasterEdges,
01751                     cutToSlaveEdges,
01752                     cutPointI,
01753                     cutPoint0,
01754                     cutPoint1
01755                 );
01756 
01757                 Pout<< "Dumping unmatched pointEdges to errorEdges.obj"
01758                     << endl;
01759 
01760                 writeOBJ
01761                 (
01762                     "errorEdges.obj",
01763                     UIndirectList<edge>
01764                     (
01765                         cutFaces().edges(),
01766                         cutFaces().pointEdges()[cutPointI]
01767                     ),
01768                     cutFaces().localPoints(),
01769                     false
01770                 );
01771 
01772                 FatalErrorIn
01773                 (
01774                     "faceCoupleInfo::subDivisionMatch"
01775                     "(const polyMesh&, const bool, const scalar)"
01776                 )   << "Problem in finding cut edges corresponding to"
01777                     << " master edge " << masterEdge
01778                     << " points:" << masterPoints[masterEdge[0]]
01779                     << ' ' << masterPoints[masterEdge[1]]
01780                     << " corresponding cut edge: (" << cutPoint0
01781                     << ' ' << cutPoint1
01782                     << abort(FatalError);
01783             }
01784 
01785             cutToMasterEdges[cutEdgeI] = masterEdgeI;
01786 
01787             cutPointI = cutEdges[cutEdgeI].otherVertex(cutPointI);
01788 
01789         } while(cutPointI != cutPoint1);
01790     }
01791 
01792     if (debug)
01793     {
01794         // Write all matched edges.
01795         writeEdges(cutToMasterEdges, cutToSlaveEdges);
01796     }
01797 
01798     // Rework cutToMasterEdges into list of points inbetween two endpoints
01799     // (cutEdgeToPoints_)
01800     setCutEdgeToPoints(cutToMasterEdges);
01801 
01802 
01803     // Now that we have marked the cut edges that lie on top of master edges
01804     // we can find cut faces that have two (or more) of these edges.
01805     // Note that there can be multiple faces having two or more matched edges
01806     // since the cut faces can form a non-manifold surface(?)
01807     // So the matching is done as an elimination process where for every
01808     // cutFace on a matched edge we store the possible master faces and
01809     // eliminate more and more until we only have one possible master face
01810     // left.
01811 
01812     if (debug)
01813     {
01814         Pout<< "subDivisionMatch :"
01815             << " matching (topological) cut faces to masterPatch" << endl;
01816     }
01817 
01818     // For every unassigned cutFaceI the possible list of master faces.
01819     Map<labelList> candidates(cutFaces().size());
01820 
01821     cutToMasterFaces_.setSize(cutFaces().size());
01822     cutToMasterFaces_ = -1;
01823 
01824     while (true)
01825     {
01826         // See if there are any edges left where there is only one remaining
01827         // candidate.
01828         label nChanged = matchEdgeFaces(cutToMasterEdges, candidates);
01829 
01830         checkMatch(cutToMasterEdges);
01831 
01832 
01833         // Now we should have found a face correspondence for every cutFace
01834         // that uses two or more cutEdges. Floodfill from these across
01835         // 'internal' cutedges (i.e. edges that do not have a master
01836         // equivalent) to determine some more cutToMasterFaces
01837         nChanged += growCutFaces(cutToMasterEdges, candidates);
01838 
01839         checkMatch(cutToMasterEdges);
01840 
01841         if (nChanged == 0)
01842         {
01843             break;
01844         }
01845     }
01846 
01847     if (debug)
01848     {
01849         Pout<< "subDivisionMatch :"
01850             << " matching (geometric) cut faces to masterPatch" << endl;
01851     }
01852 
01853     // Above should have matched all cases fully. Cannot prove this yet so
01854     // do any remaining unmatched edges by a geometric test.
01855     while (true)
01856     {
01857         // See if there are any edges left where there is only one remaining
01858         // candidate.
01859         label nChanged = geometricMatchEdgeFaces(candidates);
01860 
01861         checkMatch(cutToMasterEdges);
01862 
01863         nChanged += growCutFaces(cutToMasterEdges, candidates);
01864 
01865         checkMatch(cutToMasterEdges);
01866 
01867         if (nChanged == 0)
01868         {
01869             break;
01870         }
01871     }
01872 
01873 
01874     // All cut faces matched?
01875     forAll(cutToMasterFaces_, cutFaceI)
01876     {
01877         if (cutToMasterFaces_[cutFaceI] == -1)
01878         {
01879             const face& cutF = cutFaces()[cutFaceI];
01880 
01881             FatalErrorIn
01882             (
01883                 "faceCoupleInfo::subDivisionMatch"
01884                 "(const polyMesh&, const bool, const scalar)"
01885             )   << "Did not match all of cutFaces to a master face" << nl
01886                 << "First unmatched cut face:" << cutFaceI << " with points:"
01887                 << UIndirectList<point>(cutFaces().points(), cutF)()
01888                 << nl
01889                 << "This usually means that the slave patch is not a"
01890                 << " subdivision of the master patch"
01891                 << abort(FatalError);
01892         }
01893     }
01894 
01895     if (debug)
01896     {
01897         Pout<< "subDivisionMatch :"
01898             << " finished matching master and slave to cut" << endl;
01899     }
01900 
01901     if (debug)
01902     {
01903         writePointsFaces();
01904     }
01905 }
01906 
01907 
01908 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
01909 
01910 // Construct from mesh data
01911 Foam::faceCoupleInfo::faceCoupleInfo
01912 (
01913     const polyMesh& masterMesh,
01914     const polyMesh& slaveMesh,
01915     const scalar absTol,
01916     const bool perfectMatch
01917 )
01918 :
01919     masterPatchPtr_(NULL),
01920     slavePatchPtr_(NULL),
01921     cutPoints_(0),
01922     cutFacesPtr_(NULL),
01923     cutToMasterFaces_(0),
01924     masterToCutPoints_(0),
01925     cutToSlaveFaces_(0),
01926     slaveToCutPoints_(0),
01927     cutEdgeToPoints_(0)
01928 {
01929     // Get faces on both meshes that are aligned.
01930     // (not ordered i.e. masterToMesh[0] does
01931     // not couple to slaveToMesh[0]. This ordering is done later on)
01932     labelList masterToMesh;
01933     labelList slaveToMesh;
01934 
01935     if (perfectMatch)
01936     {
01937         // Identical faces so use tight face-centre comparison.
01938         findPerfectMatchingFaces
01939         (
01940             masterMesh,
01941             slaveMesh,
01942             absTol,
01943             masterToMesh,
01944             slaveToMesh
01945         );
01946     }
01947     else
01948     {
01949         // Slave subdivision of master so use 'nearest'. Bit dodgy, especially
01950         // with using absTol (which is quite small)
01951         findSlavesCoveringMaster
01952         (
01953             masterMesh,
01954             slaveMesh,
01955             absTol,
01956             masterToMesh,
01957             slaveToMesh
01958         );
01959     }
01960 
01961     // Construct addressing engines for both sides
01962     masterPatchPtr_.reset
01963     (
01964         new indirectPrimitivePatch
01965         (
01966             IndirectList<face>(masterMesh.faces(), masterToMesh),
01967             masterMesh.points()
01968         )
01969     );
01970 
01971     slavePatchPtr_.reset
01972     (
01973         new indirectPrimitivePatch
01974         (
01975             IndirectList<face>(slaveMesh.faces(), slaveToMesh),
01976             slaveMesh.points()
01977         )
01978     );
01979 
01980 
01981     if (perfectMatch)
01982     {
01983         // Faces are perfectly aligned but probably not ordered.
01984         perfectPointMatch(absTol, false);
01985     }
01986     else
01987     {
01988         // Slave faces are subdivision of master face. Faces not ordered.
01989         subDivisionMatch(slaveMesh, false, absTol);
01990     }
01991 
01992     if (debug)
01993     {
01994         writePointsFaces();
01995     }
01996 }
01997 
01998 
01999 // Slave is subdivision of master patch.
02000 // (so -both cover the same area -all of master points are present in slave)
02001 Foam::faceCoupleInfo::faceCoupleInfo
02002 (
02003     const polyMesh& masterMesh,
02004     const labelList& masterAddressing,
02005     const polyMesh& slaveMesh,
02006     const labelList& slaveAddressing,
02007     const scalar absTol,
02008     const bool perfectMatch,
02009     const bool orderedFaces,
02010     const bool patchDivision
02011 )
02012 :
02013     masterPatchPtr_
02014     (
02015         new indirectPrimitivePatch
02016         (
02017             IndirectList<face>(masterMesh.faces(), masterAddressing),
02018             masterMesh.points()
02019         )
02020     ),
02021     slavePatchPtr_
02022     (
02023         new indirectPrimitivePatch
02024         (
02025             IndirectList<face>(slaveMesh.faces(), slaveAddressing),
02026             slaveMesh.points()
02027         )
02028     ),
02029     cutPoints_(0),
02030     cutFacesPtr_(NULL),
02031     cutToMasterFaces_(0),
02032     masterToCutPoints_(0),
02033     cutToSlaveFaces_(0),
02034     slaveToCutPoints_(0),
02035     cutEdgeToPoints_(0)
02036 {
02037     if (perfectMatch && (masterAddressing.size() != slaveAddressing.size()))
02038     {
02039         FatalErrorIn
02040         (
02041             "faceCoupleInfo::faceCoupleInfo(const primitiveMesh&"
02042             ", const primitiveMesh&, const scalar, const bool"
02043         )   << "Perfect match specified but number of master and slave faces"
02044             << " differ." << endl
02045             << "master:" << masterAddressing.size()
02046             << "  slave:" << slaveAddressing.size()
02047             << abort(FatalError);
02048     }
02049 
02050     if
02051     (
02052         masterAddressing.size()
02053      && min(masterAddressing) < masterMesh.nInternalFaces()
02054     )
02055     {
02056         FatalErrorIn
02057         (
02058             "faceCoupleInfo::faceCoupleInfo(const primitiveMesh&"
02059             ", const primitiveMesh&, const scalar, const bool"
02060         )   << "Supplied internal face on master mesh to couple." << nl
02061             << "Faces to be coupled have to be boundary faces."
02062             << abort(FatalError);
02063     }
02064     if
02065     (
02066         slaveAddressing.size()
02067      && min(slaveAddressing) < slaveMesh.nInternalFaces()
02068     )
02069     {
02070         FatalErrorIn
02071         (
02072             "faceCoupleInfo::faceCoupleInfo(const primitiveMesh&"
02073             ", const primitiveMesh&, const scalar, const bool"
02074         )   << "Supplied internal face on slave mesh to couple." << nl
02075             << "Faces to be coupled have to be boundary faces."
02076             << abort(FatalError);
02077     }
02078 
02079 
02080     if (perfectMatch)
02081     {
02082         perfectPointMatch(absTol, orderedFaces);
02083     }
02084     else
02085     {
02086         // Slave faces are subdivision of master face. Faces not ordered.
02087         subDivisionMatch(slaveMesh, patchDivision, absTol);
02088     }
02089 
02090     if (debug)
02091     {
02092         writePointsFaces();
02093     }
02094 }
02095 
02096 
02097 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
02098 
02099 Foam::faceCoupleInfo::~faceCoupleInfo()
02100 {}
02101 
02102 
02103 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
02104 
02105 Foam::labelList Foam::faceCoupleInfo::faceLabels(const polyPatch& pp)
02106 {
02107     labelList faces(pp.size());
02108 
02109     label faceI = pp.start();
02110 
02111     forAll(pp, i)
02112     {
02113         faces[i] = faceI++;
02114     }
02115     return faces;
02116 }
02117 
02118 
02119 Foam::Map<Foam::label> Foam::faceCoupleInfo::makeMap(const labelList& lst)
02120 {
02121     Map<label> map(lst.size());
02122 
02123     forAll(lst, i)
02124     {
02125         if (lst[i] != -1)
02126         {
02127             map.insert(i, lst[i]);
02128         }
02129     }
02130     return map;
02131 }
02132 
02133 
02134 Foam::Map<Foam::labelList> Foam::faceCoupleInfo::makeMap
02135 (
02136     const labelListList& lst
02137 )
02138 {
02139     Map<labelList> map(lst.size());
02140 
02141     forAll(lst, i)
02142     {
02143         if (lst[i].size())
02144         {
02145             map.insert(i, lst[i]);
02146         }
02147     }
02148     return map;
02149 }
02150 
02151 
02152 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines