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

booleanSurface.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 Description
00025 
00026 \*---------------------------------------------------------------------------*/
00027 
00028 #include "booleanSurface.H"
00029 #include <meshTools/intersectedSurface.H>
00030 #include <meshTools/orientedSurface.H>
00031 #include <meshTools/triSurfaceSearch.H>
00032 #include <OpenFOAM/OFstream.H>
00033 #include <meshTools/treeBoundBox.H>
00034 #include <meshTools/meshTools.H>
00035 
00036 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00037 
00038 defineTypeNameAndDebug(Foam::booleanSurface, 0);
00039 
00040 
00041 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00042 
00043 // Check whether at least one of faces connected to the intersection has been
00044 // marked.
00045 void Foam::booleanSurface::checkIncluded
00046 (
00047     const intersectedSurface& surf,
00048     const labelList& faceZone,
00049     const label includedFace
00050 )
00051 {
00052     forAll(surf.intersectionEdges(), intEdgeI)
00053     {
00054         label edgeI = surf.intersectionEdges()[intEdgeI];
00055 
00056         const labelList& myFaces = surf.edgeFaces()[edgeI];
00057 
00058         bool usesIncluded = false;
00059 
00060         forAll(myFaces, myFaceI)
00061         {
00062             if (faceZone[myFaces[myFaceI]] == faceZone[includedFace])
00063             {
00064                 usesIncluded = true;
00065 
00066                 break;
00067             }
00068         }
00069 
00070         if (!usesIncluded)
00071         {
00072             FatalErrorIn
00073             (
00074                 "booleanSurface::checkIncluded(const intersectedSurface&"
00075                 ", const labelList&, const label)"
00076             )   << "None of the faces reachable from face " << includedFace
00077                 << " connects to the intersection."
00078                 << exit(FatalError);
00079         }
00080     }
00081 }
00082 
00083 
00084 // Linear lookup
00085 Foam::label Foam::booleanSurface::index
00086 (
00087     const labelList& elems,
00088     const label elem
00089 )
00090 {
00091     forAll(elems, elemI)
00092     {
00093         if (elems[elemI] == elem)
00094         {
00095             return elemI;
00096         }
00097     }
00098     return -1;
00099 }
00100 
00101 
00102 Foam::label Foam::booleanSurface::findEdge
00103 (
00104     const edgeList& edges,
00105     const labelList& edgeLabels,
00106     const edge& e
00107 )
00108 {
00109     forAll(edgeLabels, edgeLabelI)
00110     {
00111         if (edges[edgeLabels[edgeLabelI]] == e)
00112         {
00113             return edgeLabels[edgeLabelI];
00114         }
00115     }
00116     FatalErrorIn
00117     (
00118         "booleanSurface::findEdge(const edgeList&, const labelList&"
00119         ", const edge&)"
00120     )   << "Cannot find edge " << e << " in edges " << edgeLabels
00121         << abort(FatalError);
00122 
00123     return -1;
00124 }
00125 
00126 
00127 // Generate combined patchList (returned). Sets patchMap to map from surf2
00128 // region numbers into combined patchList
00129 Foam::geometricSurfacePatchList Foam::booleanSurface::mergePatches
00130 (
00131     const triSurface& surf1,
00132     const triSurface& surf2,
00133     labelList& patchMap2
00134 )
00135 {
00136     // Size too big.
00137     geometricSurfacePatchList combinedPatches
00138     (
00139         surf1.patches().size()
00140       + surf2.patches().size()
00141     );
00142 
00143     // Copy all patches of surf1
00144     label combinedPatchI = 0;
00145     forAll(surf1.patches(), patchI)
00146     {
00147         combinedPatches[combinedPatchI++] = surf1.patches()[patchI];
00148     }
00149 
00150     // (inefficiently) add unique patches from surf2
00151     patchMap2.setSize(surf2.patches().size());
00152 
00153     forAll(surf2.patches(), patch2I)
00154     {
00155         label index = -1;
00156 
00157         forAll(surf1.patches(), patch1I)
00158         {
00159             if (surf1.patches()[patch1I] == surf2.patches()[patch2I])
00160             {
00161                 index = patch1I;
00162 
00163                 break;
00164             }
00165         }
00166 
00167         if (index == -1)
00168         {
00169             combinedPatches[combinedPatchI] = surf2.patches()[patch2I];
00170             patchMap2[patch2I] = combinedPatchI;
00171             combinedPatchI++;
00172         }
00173         else
00174         {
00175             patchMap2[patch2I] = index;
00176         }
00177     }
00178 
00179     combinedPatches.setSize(combinedPatchI);
00180 
00181     return combinedPatches;
00182 }
00183 
00184 
00185 void Foam::booleanSurface::propagateEdgeSide
00186 (
00187     const triSurface& surf,
00188     const label prevVert0,
00189     const label prevFaceI,
00190     const label prevState,
00191     const label edgeI,
00192     labelList& side
00193 )
00194 {
00195     const labelList& eFaces = surf.sortedEdgeFaces()[edgeI];
00196 
00197     // Simple case. Propagate side.
00198     if (eFaces.size() == 2)
00199     {
00200         forAll(eFaces, eFaceI)
00201         {
00202             propagateSide
00203             (
00204                 surf,
00205                 prevState,
00206                 eFaces[eFaceI],
00207                 side
00208             );
00209         }
00210     }
00211 
00212 
00213     if (((eFaces.size() % 2) == 1) && (eFaces.size() != 1))
00214     {
00215         FatalErrorIn
00216         (
00217             "booleanSurface::propagateEdgeSide(const triSurface&,"
00218             "const label, const label, const label, const label,"
00219             " labelList&)"
00220         )   << "Don't know how to handle edges with odd number of faces"
00221             << endl
00222             << "edge:" << edgeI << " vertices:" << surf.edges()[edgeI]
00223             << " coming from face:" << prevFaceI
00224             << " edgeFaces:" << eFaces << abort(FatalError);
00225     }
00226 
00227 
00228     // Get position of face in edgeFaces
00229     label ind = index(eFaces, prevFaceI);
00230 
00231     // Determine orientation of faces around edge prevVert0
00232     // (might be opposite of edge)
00233     const edge& e = surf.edges()[edgeI];
00234 
00235     // Get next face to include
00236     label nextInd;
00237     label prevInd;
00238 
00239     if (e.start() == prevVert0)
00240     {
00241         // Edge (and hence eFaces) in same order as prevVert0.
00242         // Take next face from sorted list
00243         nextInd = eFaces.fcIndex(ind);
00244         prevInd = eFaces.rcIndex(ind);
00245     }
00246     else
00247     {
00248         // Take previous face from sorted neighbours
00249         nextInd = eFaces.rcIndex(ind);
00250         prevInd = eFaces.fcIndex(ind);
00251     }
00252 
00253 
00254     if (prevState == OUTSIDE)
00255     {
00256         // Coming from outside. nextInd is outside, rest is inside.
00257 
00258         forAll(eFaces, eFaceI)
00259         {
00260             if (eFaceI != ind)
00261             {
00262                 label nextState;
00263 
00264                 if (eFaceI == nextInd)
00265                 {
00266                     nextState = OUTSIDE;
00267                 }
00268                 else
00269                 {
00270                     nextState = INSIDE;
00271                 }
00272 
00273                 propagateSide
00274                 (
00275                     surf,
00276                     nextState,
00277                     eFaces[eFaceI],
00278                     side
00279                 );
00280             }
00281         }
00282     }
00283     else
00284     {
00285         // Coming from inside. prevInd is inside as well, rest is outside.
00286 
00287         forAll(eFaces, eFaceI)
00288         {
00289             if (eFaceI != ind)
00290             {
00291                 label nextState;
00292 
00293                 if (eFaceI == prevInd)
00294                 {
00295                     nextState = INSIDE;
00296                 }
00297                 else
00298                 {
00299                     nextState = OUTSIDE;
00300                 }
00301 
00302                 propagateSide
00303                 (
00304                     surf,
00305                     nextState,
00306                     eFaces[eFaceI],
00307                     side
00308                 );
00309             }
00310         }
00311     }
00312 }
00313 
00314 
00315 // Face-edge walk. Determines inside/outside for all faces connected to an edge.
00316 void Foam::booleanSurface::propagateSide
00317 (
00318     const triSurface& surf,
00319     const label prevState,
00320     const label faceI,
00321     labelList& side
00322 )
00323 {
00324     if (side[faceI] == UNVISITED)
00325     {
00326         side[faceI] = prevState;
00327 
00328         const labelledTri& tri = surf[faceI];
00329 
00330         // Get copy of face labels
00331         label a = tri[0];
00332         label b = tri[1];
00333         label c = tri[2];
00334 
00335         // Go and visit my edges' face-neighbours.
00336 
00337         const labelList& myEdges = surf.faceEdges()[faceI];
00338 
00339         label edgeAB = findEdge(surf.edges(), myEdges, edge(a, b));
00340 
00341         propagateEdgeSide
00342         (
00343             surf,
00344             a,
00345             faceI,
00346             prevState,
00347             edgeAB,
00348             side
00349         );
00350 
00351         label edgeBC = findEdge(surf.edges(), myEdges, edge(b, c));
00352 
00353         propagateEdgeSide
00354         (
00355             surf,
00356             b,
00357             faceI,
00358             prevState,
00359             edgeBC,
00360             side
00361         );
00362 
00363         label edgeCA = findEdge(surf.edges(), myEdges, edge(c, a));
00364 
00365         propagateEdgeSide
00366         (
00367             surf,
00368             c,
00369             faceI,
00370             prevState,
00371             edgeCA,
00372             side
00373         );
00374     }
00375 }
00376 
00377 
00378 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00379 
00380 // Null constructor
00381 Foam::booleanSurface::booleanSurface()
00382 :
00383     triSurface()
00384 {}
00385 
00386 
00387 // Construct from surfaces and face to include for every surface
00388 Foam::booleanSurface::booleanSurface
00389 (
00390     const triSurface& surf1,
00391     const triSurface& surf2,
00392     const surfaceIntersection& inter,
00393     const label includeFace1,
00394     const label includeFace2
00395 )
00396 :
00397     triSurface(),
00398     faceMap_()
00399 {
00400     if (debug)
00401     {
00402         Pout<< "booleanSurface : Generating intersected surface for surf1"
00403             << endl;
00404     }
00405 
00406     // Add intersection to surface1 (retriangulates cut faces)
00407     intersectedSurface cutSurf1(surf1, true, inter);
00408 
00409 
00410     if (debug)
00411     {
00412         Pout<< "booleanSurface : Generated cutSurf1: " << endl;
00413         cutSurf1.writeStats(Pout);
00414 
00415         Pout<< "Writing to file cutSurf1.ftr" << endl;
00416         OFstream cutSurf1Stream("cutSurf1.ftr");
00417         cutSurf1.write(cutSurf1Stream);
00418     }
00419 
00420     if (debug)
00421     {
00422         Pout<< "booleanSurface : Generating intersected surface for surf2"
00423             << endl;
00424     }
00425 
00426     // Add intersection to surface2
00427     intersectedSurface cutSurf2(surf2, false, inter);
00428 
00429     if (debug)
00430     {
00431         Pout<< "booleanSurface : Generated cutSurf2: " << endl;
00432         cutSurf2.writeStats(Pout);
00433 
00434         Pout<< "Writing to file cutSurf2.ftr" << endl;
00435         OFstream cutSurf2Stream("cutSurf2.ftr");
00436         cutSurf2.write(cutSurf2Stream);
00437     }
00438 
00439 
00440     // Find (first) face of cutSurf1 that originates from includeFace1
00441     label cutSurf1FaceI = index(cutSurf1.faceMap(), includeFace1);
00442 
00443     if (debug)
00444     {
00445         Pout<< "cutSurf1 : starting to fill from face:" << cutSurf1FaceI
00446             << endl;
00447     }
00448 
00449     if (cutSurf1FaceI == -1)
00450     {
00451         FatalErrorIn
00452         (
00453            "booleanSurface(const triSurfaceSearch&"
00454             ", const label, const triSurfaceSearch&, const label)"
00455         )   << "Did not find face with label " << includeFace1
00456             << " in intersectedSurface."
00457             << exit(FatalError);
00458     }
00459 
00460     // Find (first) face of cutSurf2 that originates from includeFace1
00461     label cutSurf2FaceI = index(cutSurf2.faceMap(), includeFace2);
00462 
00463     if (debug)
00464     {
00465         Pout<< "cutSurf2 : starting to fill from face:" << cutSurf2FaceI
00466             << endl;
00467     }
00468     if (cutSurf2FaceI == -1)
00469     {
00470         FatalErrorIn
00471         (
00472            "booleanSurface(const triSurfaceSearch&"
00473             ", const label, const triSurfaceSearch&, const label)"
00474         )   << "Did not find face with label " << includeFace2
00475             << " in intersectedSurface."
00476             << exit(FatalError);
00477     }
00478 
00479 
00480     //
00481     // Mark faces of cutSurf1 that need to be kept by walking from includeFace1
00482     // without crossing any edges of the intersection.
00483     //
00484 
00485     // Mark edges on intersection
00486     const labelList& int1Edges = cutSurf1.intersectionEdges();
00487 
00488     boolList isIntersectionEdge1(cutSurf1.nEdges(), false);
00489     forAll(int1Edges, intEdgeI)
00490     {
00491         label edgeI = int1Edges[intEdgeI];
00492         isIntersectionEdge1[edgeI] = true;
00493     }
00494 
00495     labelList faceZone1;
00496     (void)cutSurf1.markZones(isIntersectionEdge1, faceZone1);
00497 
00498 
00499     // Check whether at least one of sides of intersection has been marked.
00500     checkIncluded(cutSurf1, faceZone1, cutSurf1FaceI);
00501 
00502     // Subset zone which includes cutSurf2FaceI
00503     boolList includedFaces1(cutSurf1.size(), false);
00504 
00505     forAll(faceZone1, faceI)
00506     {
00507         if (faceZone1[faceI] == faceZone1[cutSurf1FaceI])
00508         {
00509             includedFaces1[faceI] = true;
00510         }
00511     }
00512 
00513     // Subset to include only interesting part
00514     labelList pointMap1;
00515     labelList faceMap1;
00516 
00517     triSurface subSurf1
00518     (
00519         cutSurf1.subsetMesh
00520         (
00521             includedFaces1,
00522             pointMap1,
00523             faceMap1
00524         )
00525     );
00526 
00527 
00528     //
00529     // Mark faces of cutSurf2 that need to be kept by walking from includeFace2
00530     // without crossing any edges of the intersection.
00531     //
00532 
00533     // Mark edges and points on intersection
00534     const labelList& int2Edges = cutSurf2.intersectionEdges();
00535 
00536     boolList isIntersectionEdge2(cutSurf2.nEdges(), false);
00537     forAll(int2Edges, intEdgeI)
00538     {
00539         label edgeI = int2Edges[intEdgeI];
00540         isIntersectionEdge2[edgeI] = true;
00541     }
00542 
00543     labelList faceZone2;
00544     (void)cutSurf2.markZones(isIntersectionEdge2, faceZone2);
00545 
00546 
00547     // Check whether at least one of sides of intersection has been marked.
00548     checkIncluded(cutSurf2, faceZone2, cutSurf2FaceI);
00549 
00550     // Subset zone which includes cutSurf2FaceI
00551     boolList includedFaces2(cutSurf2.size(), false);
00552 
00553     forAll(faceZone2, faceI)
00554     {
00555         if (faceZone2[faceI] == faceZone2[cutSurf2FaceI])
00556         {
00557             includedFaces2[faceI] = true;
00558         }
00559     }
00560 
00561     labelList pointMap2;
00562     labelList faceMap2;
00563 
00564     triSurface subSurf2
00565     (
00566         cutSurf2.subsetMesh
00567         (
00568             includedFaces2,
00569             pointMap2,
00570             faceMap2
00571         )
00572     );
00573 
00574 
00575     //
00576     // Now match up the corresponding points on the intersection. The
00577     // intersectedSurfaces will have the points resulting from the
00578     // intersection last in their points and in the same
00579     // order so we can use the pointMaps from the subsets to find them.
00580     //
00581     // We keep the vertices on the first surface and renumber those on the
00582     // second one.
00583 
00584 
00585     //
00586     // points
00587     //
00588     pointField combinedPoints
00589     (
00590         subSurf1.nPoints()
00591       + subSurf2.nPoints()
00592       - (cutSurf2.nPoints() - cutSurf2.nSurfacePoints())
00593     );
00594 
00595     // Copy points from subSurf1 and remember the labels of the ones in
00596     // the intersection
00597     labelList intersectionLabels(cutSurf1.nPoints() - cutSurf1.nSurfacePoints());
00598 
00599     label combinedPointI = 0;
00600 
00601     forAll(subSurf1.points(), pointI)
00602     {
00603         // Label in cutSurf
00604         label cutSurfPointI = pointMap1[pointI];
00605 
00606         if (!cutSurf1.isSurfacePoint(cutSurfPointI))
00607         {
00608             // Label in original intersection is equal to the cutSurfPointI
00609 
00610             // Remember label in combinedPoints for intersection point.
00611             intersectionLabels[cutSurfPointI] = combinedPointI;
00612         }
00613 
00614         // Copy point
00615         combinedPoints[combinedPointI++] = subSurf1.points()[pointI];
00616     }
00617 
00618     // Append points from subSurf2 (if they are not intersection points)
00619     // and construct mapping
00620     labelList pointMap(subSurf2.nPoints());
00621 
00622     forAll(subSurf2.points(), pointI)
00623     {
00624         // Label in cutSurf
00625         label cutSurfPointI = pointMap2[pointI];
00626 
00627         if (!cutSurf2.isSurfacePoint(cutSurfPointI))
00628         {
00629             // Lookup its label in combined point list.
00630             pointMap[pointI] = intersectionLabels[cutSurfPointI];
00631         }
00632         else
00633         {
00634             pointMap[pointI] = combinedPointI;
00635 
00636             combinedPoints[combinedPointI++] = subSurf2.points()[pointI];
00637         }
00638     }
00639 
00640 
00641     //
00642     // patches
00643     //
00644 
00645     labelList patchMap2;
00646 
00647     geometricSurfacePatchList combinedPatches
00648     (
00649         mergePatches
00650         (
00651             surf1,
00652             surf2,
00653             patchMap2
00654         )
00655     );
00656 
00657 
00658     //
00659     // faces
00660     //
00661 
00662     List<labelledTri> combinedFaces(subSurf1.size() + subSurf2.size());
00663 
00664     faceMap_.setSize(combinedFaces.size());
00665 
00666     // Copy faces from subSurf1. No need for renumbering.
00667     label combinedFaceI = 0;
00668     forAll(subSurf1, faceI)
00669     {
00670         faceMap_[combinedFaceI] = faceMap1[faceI];
00671         combinedFaces[combinedFaceI++] = subSurf1[faceI];
00672     }
00673 
00674     // Copy and renumber faces from subSurf2.
00675     forAll(subSurf2, faceI)
00676     {
00677         const labelledTri& f = subSurf2[faceI];
00678 
00679         faceMap_[combinedFaceI] = -faceMap2[faceI]-1;
00680 
00681         combinedFaces[combinedFaceI++] =
00682             labelledTri
00683             (
00684                 pointMap[f[0]],
00685                 pointMap[f[1]],
00686                 pointMap[f[2]],
00687                 patchMap2[f.region()]
00688             );
00689     }
00690 
00691     triSurface::operator=
00692     (
00693         triSurface
00694         (
00695             combinedFaces,
00696             combinedPatches,
00697             combinedPoints
00698         )
00699     );
00700 }
00701 
00702 
00703 // Construct from surfaces and boolean operation
00704 Foam::booleanSurface::booleanSurface
00705 (
00706     const triSurface& surf1,
00707     const triSurface& surf2,
00708     const surfaceIntersection& inter,
00709     const label booleanOp
00710 )
00711 :
00712     triSurface(),
00713     faceMap_()
00714 {
00715     if (debug)
00716     {
00717         Pout<< "booleanSurface : Testing surf1 and surf2" << endl;
00718 
00719         {
00720             const labelListList& edgeFaces = surf1.edgeFaces();
00721 
00722             forAll(edgeFaces, edgeI)
00723             {
00724                 const labelList& eFaces = edgeFaces[edgeI];
00725 
00726                 if (eFaces.size() == 1)
00727                 {
00728                     WarningIn("booleanSurface::booleanSurface")
00729                         << "surf1 is open surface at edge " << edgeI
00730                         << " verts:" << surf1.edges()[edgeI]
00731                         << " connected to faces " << eFaces << endl;
00732                 }
00733             }
00734         }
00735         {
00736             const labelListList& edgeFaces = surf2.edgeFaces();
00737 
00738             forAll(edgeFaces, edgeI)
00739             {
00740                 const labelList& eFaces = edgeFaces[edgeI];
00741 
00742                 if (eFaces.size() == 1)
00743                 {
00744                     WarningIn("booleanSurface::booleanSurface")
00745                         << "surf2 is open surface at edge " << edgeI
00746                         << " verts:" << surf2.edges()[edgeI]
00747                         << " connected to faces " << eFaces << endl;
00748                 }
00749             }
00750         }
00751     }
00752 
00753 
00754     //
00755     // Surface 1
00756     //
00757 
00758     if (debug)
00759     {
00760         Pout<< "booleanSurface : Generating intersected surface for surf1"
00761             << endl;
00762     }
00763 
00764     // Add intersection to surface1 (retriangulates cut faces)
00765     intersectedSurface cutSurf1(surf1, true, inter);
00766 
00767     if (debug)
00768     {
00769         Pout<< "booleanSurface : Generated cutSurf1: " << endl;
00770         cutSurf1.writeStats(Pout);
00771 
00772         Pout<< "Writing to file cutSurf1.ftr" << endl;
00773         OFstream cutSurf1Stream("cutSurf1.ftr");
00774         cutSurf1.write(cutSurf1Stream);
00775     }
00776 
00777 
00778     //
00779     // Surface 2
00780     //
00781 
00782     if (debug)
00783     {
00784         Pout<< "booleanSurface : Generating intersected surface for surf2"
00785             << endl;
00786     }
00787 
00788     // Add intersection to surface2
00789     intersectedSurface cutSurf2(surf2, false, inter);
00790 
00791     if (debug)
00792     {
00793         Pout<< "booleanSurface : Generated cutSurf2: " << endl;
00794         cutSurf2.writeStats(Pout);
00795 
00796         Pout<< "Writing to file cutSurf2.ftr" << endl;
00797         OFstream cutSurf2Stream("cutSurf2.ftr");
00798         cutSurf2.write(cutSurf2Stream);
00799     }
00800 
00801 
00802     //
00803     // patches
00804     //
00805 
00806     labelList patchMap2;
00807 
00808     geometricSurfacePatchList combinedPatches
00809     (
00810         mergePatches
00811         (
00812             surf1,
00813             surf2,
00814             patchMap2
00815         )
00816     );
00817 
00818 
00819     //
00820     // Now match up the corresponding points on the intersection. The
00821     // intersectedSurfaces will have the points resulting from the
00822     // intersection first in their points and in the same
00823     // order
00824     //
00825     // We keep the vertices on the first surface and renumber those on the
00826     // second one.
00827 
00828     pointField combinedPoints(cutSurf1.nPoints() + cutSurf2.nSurfacePoints());
00829 
00830     // Copy all points from 1 and non-intersection ones from 2.
00831 
00832     label combinedPointI = 0;
00833 
00834     forAll(cutSurf1.points(), pointI)
00835     {
00836         combinedPoints[combinedPointI++] = cutSurf1.points()[pointI];
00837     }
00838 
00839     for
00840     (
00841         label pointI = 0;
00842         pointI < cutSurf2.nSurfacePoints();
00843         pointI++
00844     )
00845     {
00846         combinedPoints[combinedPointI++] = cutSurf2.points()[pointI];
00847     }
00848 
00849     // Point order is now
00850     // - 0.. cutSurf1.nSurfacePoints : original surf1 points
00851     // -  .. cutSurf1.nPoints        : intersection points
00852     // -  .. cutSurf2.nSurfacePoints : original surf2 points
00853 
00854     if (debug)
00855     {
00856         Pout<< "booleanSurface : generated points:" << nl
00857             << "    0 .. " << cutSurf1.nSurfacePoints()-1
00858             << " : original surface1"
00859             << nl
00860             << "    " << cutSurf1.nSurfacePoints()
00861             << " .. " << cutSurf1.nPoints()-1
00862             << " : intersection points"
00863             << nl
00864             << "    " << cutSurf1.nPoints() << " .. "
00865             << cutSurf2.nSurfacePoints()-1
00866             << " : surface2 points"
00867             << nl
00868             << endl;
00869     }
00870 
00871     // Copy faces. Faces from surface 1 keep vertex numbering and region info.
00872     // Faces from 2 get vertices and region renumbered.
00873     List<labelledTri> combinedFaces(cutSurf1.size() + cutSurf2.size());
00874 
00875     label combinedFaceI = 0;
00876 
00877     forAll(cutSurf1, faceI)
00878     {
00879         combinedFaces[combinedFaceI++] = cutSurf1[faceI];
00880     }
00881 
00882     forAll(cutSurf2, faceI)
00883     {
00884         labelledTri& combinedTri = combinedFaces[combinedFaceI++];
00885 
00886         const labelledTri& tri = cutSurf2[faceI];
00887 
00888         forAll(tri, fp)
00889         {
00890             if (cutSurf2.isSurfacePoint(tri[fp]))
00891             {
00892                 // Renumber. Surface2 points are after ones from surf 1.
00893                 combinedTri[fp] = tri[fp] + cutSurf1.nPoints();
00894             }
00895             else
00896             {
00897                 // Is intersection.
00898                 combinedTri[fp] =
00899                     tri[fp]
00900                   - cutSurf2.nSurfacePoints()
00901                   + cutSurf1.nSurfacePoints();
00902             }
00903         }
00904         combinedTri.region() = patchMap2[tri.region()];
00905     }
00906 
00907 
00908     // Now we have surface in combinedFaces and combinedPoints. Use
00909     // booleanOp to determine which part of what to keep.
00910 
00911     // Construct addressing for whole part.
00912     triSurface combinedSurf
00913     (
00914         combinedFaces,
00915         combinedPatches,
00916         combinedPoints
00917     );
00918 
00919     if (debug)
00920     {
00921         Pout<< "booleanSurface : Generated combinedSurf: " << endl;
00922         combinedSurf.writeStats(Pout);
00923 
00924         Pout<< "Writing to file combinedSurf.ftr" << endl;
00925         OFstream combinedSurfStream("combinedSurf.ftr");
00926         combinedSurf.write(combinedSurfStream);
00927     }
00928 
00929 
00930     if (booleanOp == booleanSurface::ALL)
00931     {
00932         // Special case: leave surface multiply connected
00933 
00934         faceMap_.setSize(combinedSurf.size());
00935 
00936         label combinedFaceI = 0;
00937 
00938         forAll(cutSurf1, faceI)
00939         {
00940             faceMap_[combinedFaceI++] = cutSurf1.faceMap()[faceI];
00941         }
00942         forAll(cutSurf2, faceI)
00943         {
00944             faceMap_[combinedFaceI++] = -cutSurf2.faceMap()[faceI] - 1;
00945         }
00946 
00947         triSurface::operator=(combinedSurf);
00948 
00949         return;
00950     }
00951 
00952 
00953     // Get outside point.
00954     point outsidePoint = 2 * treeBoundBox(combinedSurf.localPoints()).span();
00955 
00956     //
00957     // Linear search for nearest point on surface.
00958     //
00959 
00960     const pointField& pts = combinedSurf.points();
00961 
00962     label minFaceI = -1;
00963     pointHit minHit(false, vector::zero, GREAT, true);
00964 
00965     forAll(combinedSurf, faceI)
00966     {
00967         const labelledTri& f = combinedSurf[faceI];
00968 
00969         pointHit curHit =
00970             triPointRef
00971             (
00972                 pts[f[0]],
00973                 pts[f[1]],
00974                 pts[f[2]]
00975             ).nearestPoint(outsidePoint);
00976 
00977         if (curHit.distance() < minHit.distance())
00978         {
00979             minHit = curHit;
00980 
00981             minFaceI = faceI;
00982         }
00983     }
00984 
00985     if (debug)
00986     {
00987         Pout<< "booleanSurface : found for point:" << outsidePoint
00988             << "  nearest face:" << minFaceI
00989             << "  nearest point:" << minHit.rawPoint()
00990             << endl;
00991     }
00992 
00993     // Visibility/side of face:
00994     //      UNVISITED: unvisited
00995     //      OUTSIDE: visible from outside
00996     //      INSIDE: invisible from outside
00997     labelList side(combinedSurf.size(), UNVISITED);
00998 
00999     // Walk face-edge-face and propagate inside/outside status.
01000     propagateSide(combinedSurf, OUTSIDE, minFaceI, side);
01001 
01002 
01003     // Depending on operation include certain faces.
01004     //  AND: faces on inside of 1 and of 2
01005     //  OR: faces on outside of 1 and of 2
01006     //  XOR: faces on outside of 1 and inside of 2
01007 
01008     boolList include(combinedSurf.size(), false);
01009 
01010     forAll(side, faceI)
01011     {
01012         if (side[faceI] == UNVISITED)
01013         {
01014             FatalErrorIn
01015             (
01016                 "booleanSurface::booleanSurface"
01017                 "(const triSurfaceSearch&, const triSurfaceSearch&"
01018                 ", const label booleanOp)"
01019             )   << "Face " << faceI << " has not been reached by walking from"
01020                 << " nearest point " << minHit.rawPoint()
01021                 << " nearest face " << minFaceI << exit(FatalError);
01022         }
01023         else if (side[faceI] == OUTSIDE)
01024         {
01025             if (booleanOp == booleanSurface::OR)
01026             {
01027                 include[faceI] = true;
01028             }
01029             else if (booleanOp == booleanSurface::AND)
01030             {
01031                 include[faceI] = false;
01032             }
01033             else    // xor
01034             {
01035                 include[faceI] = (faceI < cutSurf1.size()); // face from surf1
01036             }
01037         }
01038         else    // inside
01039         {
01040             if (booleanOp == booleanSurface::OR)
01041             {
01042                 include[faceI] = false;
01043             }
01044             else if (booleanOp == booleanSurface::AND)
01045             {
01046                 include[faceI] = true;
01047             }
01048             else    // xor
01049             {
01050                 include[faceI] = (faceI >= cutSurf1.size()); // face from surf2
01051             }
01052         }
01053     }
01054 
01055     // Create subsetted surface
01056     labelList subToCombinedPoint;
01057     labelList subToCombinedFace;
01058     triSurface subSurf
01059     (
01060         combinedSurf.subsetMesh
01061         (
01062             include,
01063             subToCombinedPoint,
01064             subToCombinedFace
01065         )
01066     );
01067 
01068     // Create face map
01069     faceMap_.setSize(subSurf.size());
01070 
01071     forAll(subToCombinedFace, faceI)
01072     {
01073         // Get label in combinedSurf
01074         label combinedFaceI = subToCombinedFace[faceI];
01075 
01076         // First faces in combinedSurf come from cutSurf1.
01077 
01078         if (combinedFaceI < cutSurf1.size())
01079         {
01080             label cutSurf1Face = combinedFaceI;
01081 
01082             faceMap_[faceI] = cutSurf1.faceMap()[cutSurf1Face];
01083         }
01084         else
01085         {
01086             label cutSurf2Face = combinedFaceI - cutSurf1.size();
01087 
01088             faceMap_[faceI] = - cutSurf2.faceMap()[cutSurf2Face] - 1;
01089         }
01090     }
01091 
01092     // Orient outwards
01093     orientedSurface outSurf(subSurf);
01094 
01095     // Assign.
01096     triSurface::operator=(outSurf);
01097 }
01098 
01099 
01100 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
01101 
01102 
01103 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines