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

intersectedSurface.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 "intersectedSurface.H"
00027 #include <meshTools/surfaceIntersection.H>
00028 #include <OpenFOAM/faceList.H>
00029 #include <triSurface/faceTriangulation.H>
00030 #include <meshTools/treeBoundBox.H>
00031 #include <OpenFOAM/OFstream.H>
00032 #include <OpenFOAM/error.H>
00033 #include <meshTools/meshTools.H>
00034 #include "edgeSurface.H"
00035 #include <OpenFOAM/DynamicList.H>
00036 #include <OpenFOAM/transform.H>
00037 
00038 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(Foam::intersectedSurface, 0);
00041 
00042 const Foam::label Foam::intersectedSurface::UNVISITED = 0;
00043 const Foam::label Foam::intersectedSurface::STARTTOEND = 1;
00044 const Foam::label Foam::intersectedSurface::ENDTOSTART = 2;
00045 const Foam::label Foam::intersectedSurface::BOTH = STARTTOEND | ENDTOSTART;
00046 
00047 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00048 
00049 // Write whole pointField and edges to stream
00050 void Foam::intersectedSurface::writeOBJ
00051 (
00052     const pointField& points,
00053     const edgeList& edges,
00054     Ostream& os
00055 )
00056 {
00057     forAll(points, pointI)
00058     {
00059         const point& pt = points[pointI];
00060 
00061         os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00062     }
00063     forAll(edges, edgeI)
00064     {
00065         const edge& e = edges[edgeI];
00066 
00067         os << "l " << e.start()+1 << ' ' << e.end()+1 << nl;
00068     }
00069 }
00070 
00071 
00072 // Write whole pointField and selected edges to stream
00073 void Foam::intersectedSurface::writeOBJ
00074 (
00075     const pointField& points,
00076     const edgeList& edges,
00077     const labelList& faceEdges,
00078     Ostream& os
00079 )
00080 {
00081     forAll(points, pointI)
00082     {
00083         const point& pt = points[pointI];
00084 
00085         os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00086     }
00087     forAll(faceEdges, i)
00088     {
00089         const edge& e = edges[faceEdges[i]];
00090 
00091         os << "l " << e.start()+1 << ' ' << e.end()+1 << nl;
00092     }
00093 }
00094 
00095 
00096 // write local points and edges to stream
00097 void Foam::intersectedSurface::writeLocalOBJ
00098 (
00099     const pointField& points,
00100     const edgeList& edges,
00101     const labelList& faceEdges,
00102     const fileName& fName
00103 )
00104 {
00105     OFstream os(fName);
00106 
00107     labelList pointMap(points.size(), -1);
00108 
00109     label maxVertI = 0;
00110 
00111     forAll(faceEdges, i)
00112     {
00113         const edge& e = edges[faceEdges[i]];
00114 
00115         forAll(e, i)
00116         {
00117             label pointI = e[i];
00118 
00119             if (pointMap[pointI] == -1)
00120             {
00121                 const point& pt = points[pointI];
00122 
00123                 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00124 
00125                 pointMap[pointI] = maxVertI++;
00126             }
00127         }
00128     }
00129 
00130     forAll(faceEdges, i)
00131     {
00132         const edge& e = edges[faceEdges[i]];
00133 
00134         os << "l " << pointMap[e.start()]+1 << ' ' << pointMap[e.end()]+1
00135             << nl;
00136     }
00137 }
00138 
00139 
00140 // Write whole pointField and face to stream
00141 void Foam::intersectedSurface::writeOBJ
00142 (
00143     const pointField& points,
00144     const face& f,
00145     Ostream& os
00146 )
00147 {
00148     forAll(points, pointI)
00149     {
00150         const point& pt = points[pointI];
00151 
00152         os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00153     }
00154 
00155     os << 'f';
00156 
00157     forAll(f, fp)
00158     {
00159         os << ' ' << f[fp]+1;
00160     }
00161     os << nl;
00162 }
00163 
00164 
00165 // Print current visited state.
00166 void Foam::intersectedSurface::printVisit
00167 (
00168     const edgeList& edges,
00169     const labelList& edgeLabels,
00170     const Map<label>& visited
00171 )
00172 {
00173     Pout<< "Visited:" << nl;
00174     forAll(edgeLabels, i)
00175     {
00176         label edgeI = edgeLabels[i];
00177 
00178         const edge& e = edges[edgeI];
00179 
00180         label stat = visited[edgeI];
00181 
00182         if (stat == UNVISITED)
00183         {
00184             Pout<< "    edge:" << edgeI << "  verts:" << e
00185                 << "  unvisited" << nl;
00186         }
00187         else if (stat == STARTTOEND)
00188         {
00189             Pout<< "    edge:" << edgeI << "  from " << e[0]
00190                 << " to " << e[1] << nl;
00191         }
00192         else if (stat == ENDTOSTART)
00193         {
00194             Pout<< "    edge:" << edgeI << "  from " << e[1]
00195                 << " to " << e[0] << nl;
00196         }
00197         else
00198         {
00199             Pout<< "    edge:" << edgeI << "  both " << e
00200                 << nl;
00201         }
00202     }
00203     Pout<< endl;
00204 }
00205 
00206 
00207 
00208 // Check if the two vertices that f0 and f1 share are in the same order on
00209 // both faces.
00210 bool Foam::intersectedSurface::sameEdgeOrder
00211 (
00212     const labelledTri& fA,
00213     const labelledTri& fB
00214 )
00215 {
00216     forAll(fA, fpA)
00217     {
00218         label fpB = findIndex(fB, fA[fpA]);
00219 
00220         if (fpB != -1)
00221         {
00222             // Get prev/next vertex on fA
00223             label vA1 = fA[(fpA + 1) % 3];
00224             label vAMin1 = fA[fpA ? fpA-1 : 2];
00225 
00226             // Get prev/next vertex on fB
00227             label vB1 = fB[(fpB + 1) % 3];
00228             label vBMin1 = fB[fpB ? fpB-1 : 2];
00229 
00230             if (vA1 == vB1 || vAMin1 == vBMin1)
00231             {
00232                 return true;
00233             }
00234             else if (vA1 == vBMin1 || vAMin1 == vB1)
00235             {
00236                 // shared vertices in opposite order.
00237                 return false;
00238             }
00239             else
00240             {
00241                 FatalErrorIn("intersectedSurface::sameEdgeOrder")
00242                     << "Triangle:" << fA << " and triangle:" << fB
00243                     << " share a point but not an edge"
00244                     << abort(FatalError);
00245             }
00246         }
00247     }
00248 
00249     FatalErrorIn("intersectedSurface::sameEdgeOrder")
00250         << "Triangle:" << fA << " and triangle:" << fB
00251         << " do not share an edge"
00252         << abort(FatalError);
00253 
00254     return false;
00255 }
00256 
00257 
00258 void Foam::intersectedSurface::incCount
00259 (
00260     Map<label>& visited,
00261     const label key,
00262     const label offset
00263 )
00264 {
00265     Map<label>::iterator iter = visited.find(key);
00266 
00267     if (iter == visited.end())
00268     {
00269         visited.insert(key, offset);
00270     }
00271     else
00272     {
00273         iter() += offset;
00274     }
00275 }
00276 
00277 
00278 // Calculate point to edge addressing for the face given by the edge
00279 // subset faceEdges. Constructs facePointEdges which for every point
00280 // gives a list of edge labels connected to it.
00281 Foam::Map<Foam::DynamicList<Foam::label> >
00282 Foam::intersectedSurface::calcPointEdgeAddressing
00283 (
00284     const edgeSurface& eSurf,
00285     const label faceI
00286 )
00287 {
00288     const pointField& points = eSurf.points();
00289     const edgeList& edges = eSurf.edges();
00290 
00291     const labelList& fEdges = eSurf.faceEdges()[faceI];
00292 
00293     Map<DynamicList<label> > facePointEdges(4*fEdges.size());
00294 
00295     forAll(fEdges, i)
00296     {
00297         label edgeI = fEdges[i];
00298 
00299         const edge& e = edges[edgeI];
00300 
00301         // Add e.start to point-edges
00302         Map<DynamicList<label> >::iterator iter =
00303             facePointEdges.find(e.start());
00304 
00305         if (iter == facePointEdges.end())
00306         {
00307             DynamicList<label> oneEdge;
00308             oneEdge.append(edgeI);
00309             facePointEdges.insert(e.start(), oneEdge);
00310         }
00311         else
00312         {
00313             iter().append(edgeI);
00314         }
00315 
00316         // Add e.end to point-edges
00317         Map<DynamicList<label> >::iterator iter2 =
00318             facePointEdges.find(e.end());
00319 
00320         if (iter2 == facePointEdges.end())
00321         {
00322             DynamicList<label> oneEdge;
00323             oneEdge.append(edgeI);
00324             facePointEdges.insert(e.end(), oneEdge);
00325         }
00326         else
00327         {
00328             iter2().append(edgeI);
00329         }
00330     }
00331 
00332     // Shrink it
00333     for
00334     (
00335         Map<DynamicList<label> >::iterator iter = facePointEdges.begin();
00336         iter != facePointEdges.end();
00337         ++iter
00338     )
00339     {
00340         iter().shrink();
00341 
00342         // Check on dangling points.
00343         if (iter().empty())
00344         {
00345             FatalErrorIn
00346             (
00347                 "intersectedSurface::calcPointEdgeAddressing"
00348                 "(const edgeSurface&, const label)"
00349             )   << "Point:" << iter.key() << " used by too few edges:"
00350                 << iter() << abort(FatalError);
00351         }
00352     }
00353 
00354     if (debug & 2)
00355     {
00356         // Print facePointEdges
00357         Pout<< "calcPointEdgeAddressing: face consisting of edges:" << endl;
00358         forAll(fEdges, i)
00359         {
00360             label edgeI = fEdges[i];
00361             const edge& e = edges[edgeI];
00362             Pout<< "    " << edgeI << ' ' << e << points[e.start()]
00363                 << points[e.end()] << endl;
00364         }
00365 
00366         Pout<< "    Constructed point-edge adressing:" << endl;
00367         for
00368         (
00369             Map<DynamicList<label> >::iterator iter = facePointEdges.begin();
00370             iter != facePointEdges.end();
00371             ++iter
00372         )
00373         {
00374             Pout<< "    vertex " << iter.key() << " is connected to edges "
00375                 << iter() << endl;
00376         }
00377         Pout<<endl;
00378     }
00379 
00380     return facePointEdges;
00381 }
00382 
00383 
00384 // Find next (triangle or cut) edge label coming from point prevVertI on
00385 // prevEdgeI doing a right handed walk (i.e. following right wall).
00386 // Note: normal is provided externally. Could be deducted from angle between
00387 // two triangle edges but these could be in line.
00388 Foam::label Foam::intersectedSurface::nextEdge
00389 (
00390     const edgeSurface& eSurf,
00391     const Map<label>& visited,
00392     const label faceI,
00393     const vector& n,
00394     const Map<DynamicList<label> >& facePointEdges,
00395     const label prevEdgeI,
00396     const label prevVertI
00397 )
00398 {
00399     const pointField& points = eSurf.points();
00400     const edgeList& edges = eSurf.edges();
00401     const labelList& fEdges = eSurf.faceEdges()[faceI];
00402 
00403 
00404     // Edges connected to prevVertI
00405     const DynamicList<label>& connectedEdges = facePointEdges[prevVertI];
00406 
00407     if (connectedEdges.size() <= 1)
00408     {
00409         // Problem. Point not connected.
00410         {
00411             Pout<< "Writing face:" << faceI << " to face.obj" << endl;
00412             OFstream str("face.obj");
00413             writeOBJ(points, edges, fEdges, str);
00414 
00415             Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
00416             writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
00417         }
00418 
00419         FatalErrorIn
00420         (
00421             "intersectedSurface::nextEdge(const pointField&, const edgeList&"
00422             ", const vector&, Map<DynamicList<label> >, const label"
00423             ", const label)"
00424         )   << "Problem: prevVertI:" << prevVertI << " on edge " << prevEdgeI
00425             << " has less than 2 connected edges."
00426             << " connectedEdges:" << connectedEdges << abort(FatalError);
00427 
00428         return -1;
00429     }
00430 
00431     if (connectedEdges.size() == 2)
00432     {
00433         // Simple case. Take other edge
00434         if (connectedEdges[0] == prevEdgeI)
00435         {
00436             if (debug & 2)
00437             {
00438                 Pout<< "Stepped from edge:" << edges[prevEdgeI]
00439                     << " to edge:" << edges[connectedEdges[1]]
00440                     << " chosen from candidates:" << connectedEdges << endl;
00441             }
00442             return connectedEdges[1];
00443         }
00444         else
00445         {
00446             if (debug & 2)
00447             {
00448                Pout<< "Stepped from edge:" << edges[prevEdgeI]
00449                    << " to edge:" << edges[connectedEdges[0]]
00450                    << " chosen from candidates:" << connectedEdges << endl;
00451             }
00452             return connectedEdges[0];
00453         }
00454     }
00455 
00456 
00457     // Multiple choices. Look at angle between edges.
00458 
00459     const edge& prevE = edges[prevEdgeI];
00460 
00461     // x-axis of coordinate system
00462     vector e0 = n ^ (points[prevE.otherVertex(prevVertI)] - points[prevVertI]);
00463     e0 /= mag(e0) + VSMALL;
00464 
00465     // Get y-axis of coordinate system
00466     vector e1 = n ^ e0;
00467 
00468     if (mag(mag(e1) - 1) > SMALL)
00469     {
00470         {
00471             Pout<< "Writing face:" << faceI << " to face.obj" << endl;
00472             OFstream str("face.obj");
00473             writeOBJ(points, edges, fEdges, str);
00474 
00475             Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
00476             writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
00477         }
00478 
00479         FatalErrorIn("intersectedSurface::nextEdge")
00480             << "Unnormalized normal e1:" << e1
00481             << " formed from cross product of e0:" << e0 << " n:" << n
00482             << abort(FatalError);
00483     }
00484 
00485 
00486     //
00487     // Determine maximum angle over all connected (unvisited) edges.
00488     //
00489 
00490     scalar maxAngle = -GREAT;
00491     label maxEdgeI = -1;
00492 
00493     forAll(connectedEdges, connI)
00494     {
00495         label edgeI = connectedEdges[connI];
00496 
00497         if (edgeI != prevEdgeI)
00498         {
00499             label stat = visited[edgeI];
00500 
00501             const edge& e = edges[edgeI];
00502 
00503             // Find out whether walk of edge from prevVert would be acceptible.
00504             if
00505             (
00506                 stat == UNVISITED
00507              || (
00508                     stat == STARTTOEND
00509                  && prevVertI == e[1]
00510                 )
00511              || (
00512                     stat == ENDTOSTART
00513                  && prevVertI == e[0]
00514                 )
00515             )
00516             {
00517                 // Calculate angle of edge with respect to base e0, e1
00518                 vector vec =
00519                     n ^ (points[e.otherVertex(prevVertI)] - points[prevVertI]);
00520 
00521                 vec /= mag(vec) + VSMALL;
00522 
00523                 scalar angle = pseudoAngle(e0, e1, vec);
00524 
00525                 if (angle > maxAngle)
00526                 {
00527                     maxAngle = angle;
00528                     maxEdgeI = edgeI;
00529                 }
00530             }
00531         }
00532     }
00533 
00534 
00535     if (maxEdgeI == -1)
00536     {
00537         // No unvisited edge found
00538         {
00539             Pout<< "Writing face:" << faceI << " to face.obj" << endl;
00540             OFstream str("face.obj");
00541             writeOBJ(points, edges, fEdges, str);
00542 
00543             Pout<< "Writing connectedEdges edges to faceEdges.obj" << endl;
00544             writeLocalOBJ(points, edges, connectedEdges, "faceEdges.obj");
00545         }
00546 
00547         FatalErrorIn
00548         (
00549             "intersectedSurface::nextEdge(const pointField&, const edgeList&"
00550             ", const Map<label>&, const vector&"
00551             ", const Map<DynamicList<label> >&"
00552             ", const label, const label"
00553         )   << "Trying to step from edge " << edges[prevEdgeI]
00554             << ", vertex " << prevVertI
00555             << " but cannot find 'unvisited' edges among candidates:"
00556             << connectedEdges
00557             << abort(FatalError);
00558     }
00559 
00560     if (debug & 2)
00561     {
00562         Pout<< "Stepped from edge:" << edges[prevEdgeI]
00563             << " to edge:" << maxEdgeI << " angle:" << edges[maxEdgeI]
00564             << " with angle:" << maxAngle
00565             << endl;
00566     }
00567 
00568     return maxEdgeI;
00569 }
00570 
00571 
00572 // Create (polygonal) face by walking full circle starting from startVertI on
00573 // startEdgeI.
00574 // Uses nextEdge(..) to do the walking. Returns face. Updates visited with
00575 // the labels of visited edges.
00576 Foam::face Foam::intersectedSurface::walkFace
00577 (
00578     const edgeSurface& eSurf,
00579     const label faceI,
00580     const vector& n,
00581     const Map<DynamicList<label> >& facePointEdges,
00582 
00583     const label startEdgeI,
00584     const label startVertI,
00585 
00586     Map<label>& visited
00587 )
00588 {
00589     const pointField& points = eSurf.points();
00590     const edgeList& edges = eSurf.edges();
00591 
00592     // Overestimate size of face
00593     face f(eSurf.faceEdges()[faceI].size());
00594 
00595     label fp = 0;
00596 
00597     label vertI = startVertI;
00598     label edgeI = startEdgeI;
00599 
00600     while(true)
00601     {
00602         const edge& e = edges[edgeI];
00603 
00604         if (debug & 2)
00605         {
00606             Pout<< "Now at:" << endl
00607                 << "    edge:" << edgeI << " vertices:" << e
00608                 << " positions:" << points[e.start()] << ' ' << points[e.end()]
00609                 << "    vertex:" << vertI << endl;
00610         }
00611 
00612         // Mark edge as visited
00613         if (e[0] == vertI)
00614         {
00615             visited[edgeI] |= STARTTOEND;
00616         }
00617         else
00618         {
00619             visited[edgeI] |= ENDTOSTART;
00620         }
00621 
00622 
00623         // Store face vertex
00624         f[fp++] = vertI;
00625 
00626         // step to other vertex
00627         vertI = e.otherVertex(vertI);
00628 
00629         if (vertI == startVertI)
00630         {
00631             break;
00632         }
00633 
00634         // step from vertex to next edge
00635         edgeI = nextEdge
00636         (
00637             eSurf,
00638             visited,
00639             faceI,
00640             n,
00641             facePointEdges,
00642             edgeI,
00643             vertI
00644         );
00645     }
00646 
00647     f.setSize(fp);
00648 
00649     return f;
00650 }
00651 
00652 
00653 void Foam::intersectedSurface::findNearestVisited
00654 (
00655     const edgeSurface& eSurf,
00656     const label faceI,
00657     const Map<DynamicList<label> >& facePointEdges,
00658     const Map<label>& pointVisited,
00659     const point& pt,
00660     const label excludePointI,
00661 
00662     label& minVertI,
00663     scalar& minDist
00664 )
00665 {
00666     minVertI = -1;
00667     minDist = GREAT;
00668 
00669     forAllConstIter(Map<label>, pointVisited, iter)
00670     {
00671         label pointI = iter.key();
00672 
00673         if (pointI != excludePointI)
00674         {
00675             label nVisits = iter();
00676 
00677             if (nVisits == 2*facePointEdges[pointI].size())
00678             {
00679                 // Fully visited (i.e. both sides of all edges)
00680                 scalar dist = mag(eSurf.points()[pointI] - pt);
00681 
00682                 if (dist < minDist)
00683                 {
00684                     minDist = dist;
00685                     minVertI = pointI;
00686                 }
00687             }
00688         }
00689     }
00690 
00691     if (minVertI == -1)
00692     {
00693         const labelList& fEdges = eSurf.faceEdges()[faceI];
00694 
00695         SeriousErrorIn("intersectedSurface::findNearestVisited")
00696             << "Dumping face edges to faceEdges.obj" << endl;
00697 
00698         writeLocalOBJ(eSurf.points(), eSurf.edges(), fEdges, "faceEdges.obj");
00699 
00700         FatalErrorIn("intersectedSurface::findNearestVisited")
00701             << "No fully visited edge found for pt " << pt
00702             << abort(FatalError);
00703     }
00704 }
00705 
00706 
00707 // Sometimes there are bunches of edges that are not connected to the
00708 // other edges in the face. This routine tries to connect the loose
00709 // edges up to the 'proper' edges by adding two extra edges between a
00710 // properly connected edge and an unconnected one. Since at this level the
00711 // adressing is purely in form of points and a cloud of edges this can
00712 // be easily done.
00713 // (edges are to existing points. Don't want to introduce new vertices here
00714 // since then also neighbouring face would have to be split)
00715 Foam::faceList Foam::intersectedSurface::resplitFace
00716 (
00717     const triSurface& surf,
00718     const label faceI,
00719     const Map<DynamicList<label> >& facePointEdges,
00720     const Map<label>& visited,
00721     edgeSurface& eSurf
00722 )
00723 {
00724     {
00725         // Dump face for debugging.
00726         Pout<< "Writing face:" << faceI << " to face.obj" << endl;
00727         OFstream str("face.obj");
00728         writeOBJ(eSurf.points(), eSurf.edges(), eSurf.faceEdges()[faceI], str);
00729     }
00730 
00731 
00732     // Count the number of times point has been visited so we
00733     // can compare number to facePointEdges.
00734     Map<label> pointVisited(2*facePointEdges.size());
00735 
00736 
00737     {
00738         OFstream str0("visitedNone.obj");
00739         OFstream str1("visitedOnce.obj");
00740         OFstream str2("visitedTwice.obj");
00741         forAll(eSurf.points(), pointI)
00742         {
00743             const point& pt = eSurf.points()[pointI];
00744 
00745             str0 << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00746             str1 << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00747             str2 << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
00748         }
00749 
00750 
00751     forAllConstIter(Map<label>, visited, iter)
00752     {
00753         label edgeI = iter.key();
00754 
00755         const edge& e = eSurf.edges()[edgeI];
00756 
00757         label stat = iter();
00758 
00759         if (stat == STARTTOEND || stat == ENDTOSTART)
00760         {
00761             incCount(pointVisited, e[0], 1);
00762             incCount(pointVisited, e[1], 1);
00763 
00764             str1 << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00765         }
00766         else if (stat == BOTH)
00767         {
00768             incCount(pointVisited, e[0], 2);
00769             incCount(pointVisited, e[1], 2);
00770 
00771             str2 << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00772         }
00773         else if (stat == UNVISITED)
00774         {
00775             incCount(pointVisited, e[0], 0);
00776             incCount(pointVisited, e[1], 0);
00777 
00778             str0 << "l " << e[0]+1 << ' ' << e[1]+1 << nl;
00779         }
00780     }
00781     }
00782 
00783 
00784     {
00785         forAllConstIter(Map<label>, pointVisited, iter)
00786         {
00787             label pointI = iter.key();
00788 
00789             label nVisits = iter();
00790 
00791             Pout<< "point:" << pointI << "  nVisited:" << nVisits
00792                 << "  pointEdges:" << facePointEdges[pointI].size() << endl;
00793         }
00794     }
00795 
00796 
00797     // Find nearest point pair where one is not fully visited and
00798     // the other is.
00799 
00800     label visitedVert0 = -1;
00801     label unvisitedVert0 = -1;
00802 
00803     {
00804         scalar minDist = GREAT;
00805 
00806         forAllConstIter(Map<DynamicList<label> >, facePointEdges, iter)
00807         {
00808             label pointI = iter.key();
00809 
00810             label nVisits = pointVisited[pointI];
00811 
00812             const DynamicList<label>& pEdges = iter();
00813 
00814             if (nVisits < 2*pEdges.size())
00815             {
00816                 // Not fully visited. Find nearest fully visited.
00817 
00818                 scalar nearDist;
00819                 label nearVertI;
00820 
00821                 findNearestVisited
00822                 (
00823                     eSurf,
00824                     faceI,
00825                     facePointEdges,
00826                     pointVisited,
00827                     eSurf.points()[pointI],
00828                     -1,                         // Do not exclude vertex
00829                     nearVertI,
00830                     nearDist
00831                 );
00832 
00833 
00834                 if (nearDist < minDist)
00835                 {
00836                     minDist = nearDist;
00837                     visitedVert0 = nearVertI;
00838                     unvisitedVert0 = pointI;
00839                 }
00840             }
00841         }
00842     }
00843 
00844 
00845     // Find second intersection.
00846     label visitedVert1 = -1;
00847     label unvisitedVert1 = -1;
00848     {
00849         scalar minDist = GREAT;
00850 
00851         forAllConstIter(Map<DynamicList<label> >, facePointEdges, iter)
00852         {
00853             label pointI = iter.key();
00854 
00855             if (pointI != unvisitedVert0)
00856             {
00857                 label nVisits = pointVisited[pointI];
00858 
00859                 const DynamicList<label>& pEdges = iter();
00860 
00861                 if (nVisits < 2*pEdges.size())
00862                 {
00863                     // Not fully visited. Find nearest fully visited.
00864 
00865                     scalar nearDist;
00866                     label nearVertI;
00867 
00868                     findNearestVisited
00869                     (
00870                         eSurf,
00871                         faceI,
00872                         facePointEdges,
00873                         pointVisited,
00874                         eSurf.points()[pointI],
00875                         visitedVert0,           // vertex to exclude
00876                         nearVertI,
00877                         nearDist
00878                     );
00879 
00880 
00881                     if (nearDist < minDist)
00882                     {
00883                         minDist = nearDist;
00884                         visitedVert1 = nearVertI;
00885                         unvisitedVert1 = pointI;
00886                     }
00887                 }
00888             }
00889         }
00890     }
00891 
00892 
00893     Pout<< "resplitFace : adding intersection from " << visitedVert0
00894         << " to " << unvisitedVert0 << endl
00895         << " and from " << visitedVert1
00896         << " to " << unvisitedVert1 << endl;
00897 
00898 
00899 //    // Add the new intersection edges to the edgeSurface
00900 //    edgeList additionalEdges(2);
00901 //    additionalEdges[0] = edge(visitedVert0, unvisitedVert0);
00902 //    additionalEdges[1] = edge(visitedVert1, unvisitedVert1);
00903 
00904     // Add the new intersection edges to the edgeSurface
00905     edgeList additionalEdges(1);
00906     additionalEdges[0] = edge(visitedVert0, unvisitedVert0);
00907 
00908     eSurf.addIntersectionEdges(faceI, additionalEdges);
00909 
00910     fileName newFName("face_" + Foam::name(faceI) + "_newEdges.obj");
00911     Pout<< "Dumping face:" << faceI << " to " << newFName << endl;
00912     writeLocalOBJ
00913     (
00914         eSurf.points(),
00915         eSurf.edges(),
00916         eSurf.faceEdges()[faceI],
00917         newFName
00918     );
00919 
00920     // Retry splitFace. Use recursion since is rare situation.
00921     return splitFace(surf, faceI, eSurf);
00922 }
00923 
00924 
00925 Foam::faceList Foam::intersectedSurface::splitFace
00926 (
00927     const triSurface& surf,
00928     const label faceI,
00929     edgeSurface& eSurf
00930 )
00931 {
00932     // Alias
00933     const pointField& points = eSurf.points();
00934     const edgeList& edges = eSurf.edges();
00935     const labelList& fEdges = eSurf.faceEdges()[faceI];
00936 
00937     // Create local (for the face only) point-edge connectivity.
00938     Map<DynamicList<label> > facePointEdges
00939     (
00940         calcPointEdgeAddressing
00941         (
00942             eSurf,
00943             faceI
00944         )
00945     );
00946 
00947     // Order in which edges have been walked. Initialize outside edges.
00948     Map<label> visited(fEdges.size()*2);
00949 
00950     forAll(fEdges, i)
00951     {
00952         label edgeI = fEdges[i];
00953 
00954         if (eSurf.isSurfaceEdge(edgeI))
00955         {
00956             // Edge is edge from original surface so an outside edge for
00957             // the current face.
00958             label surfEdgeI = eSurf.parentEdge(edgeI);
00959 
00960             label owner = surf.edgeOwner()[surfEdgeI];
00961 
00962             if
00963             (
00964                 owner == faceI
00965              || sameEdgeOrder
00966                 (
00967                     surf.localFaces()[owner],
00968                     surf.localFaces()[faceI]
00969                 )
00970             )
00971             {
00972                 // Edge is in same order as current face.
00973                 // Mark off the opposite order.
00974                 visited.insert(edgeI, ENDTOSTART);
00975             }
00976             else
00977             {
00978                 // Edge is in same order as current face.
00979                 // Mark off the opposite order.
00980                 visited.insert(edgeI, STARTTOEND);
00981             }
00982         }
00983         else
00984         {
00985             visited.insert(edgeI, UNVISITED);
00986         }
00987     }
00988 
00989 
00990 
00991     // Storage for faces
00992     DynamicList<face> faces(fEdges.size());
00993 
00994     while (true)
00995     {
00996         // Find starting edge:
00997         // - unvisted triangle edge
00998         // - once visited intersection edge
00999         // Give priority to triangle edges.
01000         label startEdgeI = -1;
01001         label startVertI = -1;
01002 
01003         forAll(fEdges, i)
01004         {
01005             label edgeI = fEdges[i];
01006 
01007             const edge& e = edges[edgeI];
01008 
01009             label stat = visited[edgeI];
01010 
01011             if (stat == STARTTOEND)
01012             {
01013                 startEdgeI = edgeI;
01014                 startVertI = e[1];
01015 
01016                 if (eSurf.isSurfaceEdge(edgeI))
01017                 {
01018                     break;
01019                 }
01020             }
01021             else if (stat == ENDTOSTART)
01022             {
01023                 startEdgeI = edgeI;
01024                 startVertI = e[0];
01025 
01026                 if (eSurf.isSurfaceEdge(edgeI))
01027                 {
01028                     break;
01029                 }
01030             }
01031         }
01032 
01033         if (startEdgeI == -1)
01034         {
01035             break;
01036         }
01037 
01038         //Pout<< "splitFace: starting walk from edge:" << startEdgeI
01039         //    << ' ' << edges[startEdgeI] << " vertex:" << startVertI << endl;
01040 
01042         //printVisit(eSurf.edges(), fEdges, visited);
01043 
01044         //{
01045         //    Pout<< "Writing face:" << faceI << " to face.obj" << endl;
01046         //    OFstream str("face.obj");
01047         //    writeOBJ(eSurf.points(), eSurf.edges(), fEdges, str);
01048         //}
01049 
01050         faces.append
01051         (
01052             walkFace
01053             (
01054                 eSurf,
01055                 faceI,
01056                 surf.faceNormals()[faceI],
01057                 facePointEdges,
01058 
01059                 startEdgeI,
01060                 startVertI,
01061 
01062                 visited
01063             )
01064         );
01065     }
01066 
01067 
01068     // Check if any unvisited edges left.
01069     forAll(fEdges, i)
01070     {
01071         label edgeI = fEdges[i];
01072 
01073         label stat = visited[edgeI];
01074 
01075         if (eSurf.isSurfaceEdge(edgeI) && stat != BOTH)
01076         {
01077             SeriousErrorIn("Foam::intersectedSurface::splitFace")
01078                 << "Dumping face edges to faceEdges.obj" << endl;
01079 
01080             writeLocalOBJ(points, edges, fEdges, "faceEdges.obj");
01081 
01082             FatalErrorIn("intersectedSurface::splitFace")
01083                << "Problem: edge " << edgeI << " vertices "
01084                 << edges[edgeI] << " on face " << faceI
01085                 << " has visited status " << stat << " from a "
01086                  << "righthanded walk along all"
01087                 << " of the triangle edges. Are the original surfaces"
01088                 << " closed and non-intersecting?"
01089                 << abort(FatalError);
01090         }
01091         else if (stat != BOTH)
01092         {
01093             //{
01094             //    Pout<< "Dumping faces so far to faces.obj" << nl
01095             //        << faces << endl;
01096             //
01097             //    OFstream str("faces.obj");
01098             //
01099             //    forAll(faces, i)
01100             //    {
01101             //        writeOBJ(points, faces[i], str);
01102             //    }
01103             //}
01104 
01105             Pout<< "** Resplitting **" << endl;
01106 
01107             // Redo face after introducing extra edge. Edge introduced
01108             // should be one nearest to any fully visited edge.
01109             return resplitFace
01110             (
01111                 surf,
01112                 faceI,
01113                 facePointEdges,
01114                 visited,
01115                 eSurf
01116             );
01117         }
01118     }
01119 
01120 
01121     // See if normal needs flipping.
01122     faces.shrink();
01123 
01124     vector n = faces[0].normal(eSurf.points());
01125 
01126     if ((n & surf.faceNormals()[faceI]) < 0)
01127     {
01128         forAll(faces, i)
01129         {
01130             reverse(faces[i]);
01131         }
01132     }
01133 
01134     return faces;
01135 }
01136 
01137 
01138 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
01139 
01140 // Null constructor
01141 Foam::intersectedSurface::intersectedSurface()
01142 :
01143     triSurface(),
01144     intersectionEdges_(0),
01145     faceMap_(0),
01146     nSurfacePoints_(0)
01147 {}
01148 
01149 
01150 // Construct from components
01151 Foam::intersectedSurface::intersectedSurface(const triSurface& surf)
01152 :
01153     triSurface(surf),
01154     intersectionEdges_(0),
01155     faceMap_(0),
01156     nSurfacePoints_(0)
01157 {}
01158 
01159 
01160 // Construct from surface and intersection
01161 Foam::intersectedSurface::intersectedSurface
01162 (
01163     const triSurface& surf,
01164     const bool isFirstSurface,
01165     const surfaceIntersection& inter
01166 )
01167 :
01168     triSurface(),
01169     intersectionEdges_(0),
01170     faceMap_(0),
01171     nSurfacePoints_(surf.nPoints())
01172 {
01173     if (inter.cutPoints().empty() && inter.cutEdges().empty())
01174     {
01175         // No intersection. Make straight copy.
01176         triSurface::operator=(surf);
01177 
01178         // Identity for face map
01179         faceMap_.setSize(size());
01180 
01181         forAll(faceMap_, faceI)
01182         {
01183             faceMap_[faceI] = faceI;
01184         }
01185         return;
01186     }
01187 
01188 
01189     // Calculate face-edge addressing for surface + intersection.
01190     edgeSurface eSurf(surf, isFirstSurface, inter);
01191 
01192     // Update point information for any removed points. (when are they removed?
01193     // - but better make sure)
01194     nSurfacePoints_ = eSurf.nSurfacePoints();
01195 
01196     // Now we have full points, edges and edge addressing for surf. Start
01197     // extracting faces and triangulate them.
01198 
01199     // Storage for new triangles of surface.
01200     DynamicList<labelledTri> newTris(eSurf.edges().size()/2);
01201 
01202     // Start in newTris for decomposed face.
01203     labelList startTriI(surf.size(), 0);
01204 
01205     forAll(surf, faceI)
01206     {
01207         startTriI[faceI] = newTris.size();
01208 
01209         if (eSurf.faceEdges()[faceI].size() != surf.faceEdges()[faceI].size())
01210         {
01211             // Face has been cut by intersection.
01212             // Cut face into multiple subfaces. Use faceEdge information
01213             // from edgeSurface only. (original triSurface 'surf' is used for
01214             // faceNormals and regionnumber only)
01215             faceList newFaces
01216             (
01217                 splitFace
01218                 (
01219                     surf,
01220                     faceI,              // current triangle
01221                     eSurf               // face-edge description of surface
01222                                         // + intersection
01223                 )
01224             );
01225             forAll(newFaces, newFaceI)
01226             {
01227                 const face& newF = newFaces[newFaceI];
01228 
01229 //                {
01230 //                    fileName fName
01231 //                    (
01232 //                        "face_"
01233 //                      + Foam::name(faceI)
01234 //                      + "_subFace_"
01235 //                      + Foam::name(newFaceI)
01236 //                      + ".obj"
01237 //                    );
01238 //                    Pout<< "Writing original face:" << faceI << " subFace:"
01239 //                        << newFaceI << " to " << fName << endl;
01240 //
01241 //                    OFstream str(fName);
01242 //
01243 //                    forAll(newF, fp)
01244 //                    {
01245 //                        meshTools::writeOBJ(str, eSurf.points()[newF[fp]]);
01246 //                    }
01247 //                    str << 'l';
01248 //                    forAll(newF, fp)
01249 //                    {
01250 //                        str << ' ' << fp+1;
01251 //                    }
01252 //                    str<< " 1" << nl;
01253 //                }
01254 
01255 
01256                 const vector& n = surf.faceNormals()[faceI];
01257                 const label region = surf[faceI].region();
01258 
01259                 faceTriangulation tris(eSurf.points(), newF, n);
01260 
01261                 forAll(tris, triI)
01262                 {
01263                     const triFace& t = tris[triI];
01264 
01265                     forAll(t, i)
01266                     {
01267                         if (t[i] < 0 || t[i] >= eSurf.points().size())
01268                         {
01269                             FatalErrorIn
01270                             (
01271                                 "intersectedSurface::intersectedSurface"
01272                             )   << "Face triangulation of face " << faceI
01273                                 << " uses points outside range 0.."
01274                                 << eSurf.points().size()-1 << endl
01275                                 << "Triangulation:"
01276                                 << tris << abort(FatalError);
01277                         }
01278                     }
01279 
01280                     newTris.append(labelledTri(t[0], t[1], t[2], region));
01281                 }
01282             }
01283         }
01284         else
01285         {
01286             // Face has not been cut at all. No need to renumber vertices since
01287             // eSurf keeps surface vertices first.
01288             newTris.append(surf.localFaces()[faceI]);
01289         }
01290     }
01291 
01292     newTris.shrink();
01293 
01294 
01295     // Construct triSurface. Note that addressing will be compact since
01296     // edgeSurface is compact.
01297     triSurface::operator=
01298     (
01299         triSurface
01300         (
01301             newTris,
01302             surf.patches(),
01303             eSurf.points()
01304         )
01305     );
01306 
01307     // Construct mapping back into original surface
01308     faceMap_.setSize(size());
01309 
01310     for (label faceI = 0; faceI < surf.size()-1; faceI++)
01311     {
01312         for (label triI = startTriI[faceI]; triI < startTriI[faceI+1]; triI++)
01313         {
01314             faceMap_[triI] = faceI;
01315         }
01316     }
01317     for (label triI = startTriI[surf.size()-1]; triI < size(); triI++)
01318     {
01319         faceMap_[triI] = surf.size()-1;
01320     }
01321 
01322 
01323     // Find edges on *this which originate from 'cuts'. (i.e. newEdgeI >=
01324     // nSurfaceEdges) Renumber edges into local triSurface numbering.
01325 
01326     intersectionEdges_.setSize(eSurf.edges().size() - eSurf.nSurfaceEdges());
01327 
01328     label intersectionEdgeI = 0;
01329 
01330     for
01331     (
01332         label edgeI = eSurf.nSurfaceEdges();
01333         edgeI < eSurf.edges().size();
01334         edgeI++
01335     )
01336     {
01337         // Get edge vertices in triSurface local numbering
01338         const edge& e = eSurf.edges()[edgeI];
01339         label surfStartI = meshPointMap()[e.start()];
01340         label surfEndI = meshPointMap()[e.end()];
01341 
01342         // Find edge connected to surfStartI which also uses surfEndI.
01343         label surfEdgeI = -1;
01344 
01345         const labelList& pEdges = pointEdges()[surfStartI];
01346 
01347         forAll(pEdges, i)
01348         {
01349             const edge& surfE = edges()[pEdges[i]];
01350 
01351             // Edge already connected to surfStart for sure. See if also
01352             // connects to surfEnd
01353             if (surfE.start() == surfEndI || surfE.end() == surfEndI)
01354             {
01355                 surfEdgeI = pEdges[i];
01356 
01357                 break;
01358             }
01359         }
01360 
01361         if (surfEdgeI != -1)
01362         {
01363             intersectionEdges_[intersectionEdgeI++] = surfEdgeI;
01364         }
01365         else
01366         {
01367             FatalErrorIn("intersectedSurface::intersectedSurface")
01368                 << "Cannot find edge among candidates " << pEdges
01369                 << " which uses points " << surfStartI
01370                 << " and " << surfEndI
01371                 << abort(FatalError);
01372         }
01373     }
01374 }
01375 
01376 
01377 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines