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

regionSide.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 "regionSide.H"
00027 #include <meshTools/meshTools.H>
00028 #include <OpenFOAM/primitiveMesh.H>
00029 
00030 
00031 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 
00036 defineTypeNameAndDebug(regionSide, 0);
00037 
00038 }
00039 
00040 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00041 
00042 // Step across edge onto other face on cell
00043 Foam::label Foam::regionSide::otherFace
00044 (
00045     const primitiveMesh& mesh,
00046     const label cellI,
00047     const label faceI,
00048     const label edgeI
00049 )
00050 {
00051     label f0I, f1I;
00052 
00053     meshTools::getEdgeFaces(mesh, cellI, edgeI, f0I, f1I);
00054 
00055     if (f0I == faceI)
00056     {
00057         return f1I;
00058     }
00059     else
00060     {
00061         return f0I;
00062     }
00063 }
00064 
00065 
00066 // Step across point to other edge on face
00067 Foam::label Foam::regionSide::otherEdge
00068 (
00069     const primitiveMesh& mesh,
00070     const label faceI,
00071     const label edgeI,
00072     const label pointI
00073 )
00074 {
00075     const edge& e = mesh.edges()[edgeI];
00076 
00077     // Get other point on edge.
00078     label freePointI = e.otherVertex(pointI);
00079 
00080     const labelList& fEdges = mesh.faceEdges()[faceI];
00081 
00082     forAll(fEdges, fEdgeI)
00083     {
00084         label otherEdgeI = fEdges[fEdgeI];
00085 
00086         const edge& otherE = mesh.edges()[otherEdgeI];
00087 
00088         if
00089         (
00090             (
00091                 otherE.start() == pointI
00092              && otherE.end() != freePointI
00093             )
00094          || (
00095                 otherE.end() == pointI
00096              && otherE.start() != freePointI
00097             )
00098         )
00099         {
00100             // otherE shares one (but not two) points with e.
00101             return otherEdgeI;
00102         }
00103     }
00104 
00105     FatalErrorIn
00106     (
00107         "regionSide::otherEdge(const primitiveMesh&, const label, const label"
00108         ", const label)"
00109     )   << "Cannot find other edge on face " << faceI << " that uses point "
00110         << pointI << " but not point " << freePointI << endl
00111         << "Edges on face:" << fEdges
00112         << " verts:" << UIndirectList<edge>(mesh.edges(), fEdges)()
00113         << " Vertices on face:"
00114         << mesh.faces()[faceI]
00115         << " Vertices on original edge:" << e << abort(FatalError);
00116 
00117     return -1;
00118 }
00119 
00120 
00121 // Step from faceI (on side cellI) to connected face & cell without crossing
00122 // fenceEdges.
00123 void Foam::regionSide::visitConnectedFaces
00124 (
00125     const primitiveMesh& mesh,
00126     const labelHashSet& region,
00127     const labelHashSet& fenceEdges,
00128     const label cellI,
00129     const label faceI,
00130     labelHashSet& visitedFace
00131 )
00132 {
00133     if (!visitedFace.found(faceI))
00134     {
00135         if (debug)
00136         {
00137             Info<< "visitConnectedFaces : cellI:" << cellI << " faceI:"
00138                 << faceI << "  isOwner:" << (cellI == mesh.faceOwner()[faceI])
00139                 << endl;
00140         }
00141 
00142         // Mark as visited
00143         visitedFace.insert(faceI);
00144 
00145         // Mark which side of face was visited.
00146         if (cellI == mesh.faceOwner()[faceI])
00147         {
00148             sideOwner_.insert(faceI);
00149         }
00150 
00151 
00152         // Visit all neighbouring faces on faceSet. Stay on this 'side' of
00153         // face by doing edge-face-cell walk.
00154         const labelList& fEdges = mesh.faceEdges()[faceI];
00155 
00156         forAll(fEdges, fEdgeI)
00157         {
00158             label edgeI = fEdges[fEdgeI];
00159 
00160             if (!fenceEdges.found(edgeI))
00161             {
00162                 // Step along faces on edge from cell to cell until
00163                 // we hit face on faceSet.
00164 
00165                 // Find face reachable from edge
00166                 label otherFaceI = otherFace(mesh, cellI, faceI, edgeI);
00167 
00168                 if (mesh.isInternalFace(otherFaceI))
00169                 {
00170                     label otherCellI = cellI;
00171 
00172                     // Keep on crossing faces/cells until back on face on
00173                     // surface
00174                     while (!region.found(otherFaceI))
00175                     {
00176                         visitedFace.insert(otherFaceI);
00177 
00178                         if (debug)
00179                         {
00180                             Info<< "visitConnectedFaces : cellI:" << cellI
00181                                 << " found insideEdgeFace:" << otherFaceI
00182                                 << endl;
00183                         }
00184 
00185 
00186                         // Cross otherFaceI into neighbouring cell
00187                         otherCellI =
00188                             meshTools::otherCell
00189                             (
00190                                 mesh,
00191                                 otherCellI,
00192                                 otherFaceI
00193                             );
00194 
00195                         otherFaceI =
00196                                 otherFace
00197                                 (
00198                                     mesh,
00199                                     otherCellI,
00200                                     otherFaceI,
00201                                     edgeI
00202                                 );
00203                     }
00204 
00205                     visitConnectedFaces
00206                     (
00207                         mesh,
00208                         region,
00209                         fenceEdges,
00210                         otherCellI,
00211                         otherFaceI,
00212                         visitedFace
00213                     );
00214                 }
00215             }
00216         }
00217     }
00218 }
00219 
00220 
00221 // From edge on face connected to point on region (regionPointI) cross
00222 // to all other edges using this point by walking across faces
00223 // Does not cross regionEdges so stays on one side
00224 // of region
00225 void Foam::regionSide::walkPointConnectedFaces
00226 (
00227     const primitiveMesh& mesh,
00228     const labelHashSet& regionEdges,
00229     const label regionPointI,
00230     const label startFaceI,
00231     const label startEdgeI,
00232     labelHashSet& visitedEdges
00233 )
00234 {
00235     // Mark as visited
00236     insidePointFaces_.insert(startFaceI);
00237 
00238     if (debug)
00239     {
00240         Info<< "walkPointConnectedFaces : regionPointI:" << regionPointI
00241             << " faceI:" << startFaceI
00242             << " edgeI:" << startEdgeI << " verts:"
00243             << mesh.edges()[startEdgeI]
00244             << endl;
00245     }
00246 
00247     // Cross faceI i.e. get edge not startEdgeI which uses regionPointI
00248     label edgeI = otherEdge(mesh, startFaceI, startEdgeI, regionPointI);
00249 
00250     if (!regionEdges.found(edgeI))
00251     {
00252         if (!visitedEdges.found(edgeI))
00253         {
00254             visitedEdges.insert(edgeI);
00255 
00256             if (debug)
00257             {
00258                 Info<< "Crossed face from "
00259                     << " edgeI:" << startEdgeI << " verts:"
00260                     << mesh.edges()[startEdgeI]
00261                     << " to edge:" << edgeI << " verts:"
00262                     << mesh.edges()[edgeI]
00263                     << endl;
00264             }
00265 
00266             // Cross edge to all faces connected to it.
00267 
00268             const labelList& eFaces = mesh.edgeFaces()[edgeI];
00269 
00270             forAll(eFaces, eFaceI)
00271             {
00272                 label faceI = eFaces[eFaceI];
00273 
00274                 walkPointConnectedFaces
00275                 (
00276                     mesh,
00277                     regionEdges,
00278                     regionPointI,
00279                     faceI,
00280                     edgeI,
00281                     visitedEdges
00282                 );
00283             }
00284         }
00285     }
00286 }
00287 
00288 
00289 // Find all faces reachable from all non-fence points and staying on
00290 // regionFaces side.
00291 void Foam::regionSide::walkAllPointConnectedFaces
00292 (
00293     const primitiveMesh& mesh,
00294     const labelHashSet& regionFaces,
00295     const labelHashSet& fencePoints
00296 )
00297 {
00298     //
00299     // Get all (internal and external) edges on region.
00300     //
00301     labelHashSet regionEdges(4*regionFaces.size());
00302 
00303     for
00304     (
00305         labelHashSet::const_iterator iter = regionFaces.begin();
00306         iter != regionFaces.end();
00307         ++iter
00308     )
00309     {
00310         label faceI = iter.key();
00311 
00312         const labelList& fEdges = mesh.faceEdges()[faceI];
00313 
00314         forAll(fEdges, fEdgeI)
00315         {
00316             regionEdges.insert(fEdges[fEdgeI]);
00317         }
00318     }
00319 
00320 
00321     //
00322     // Visit all internal points on surface.
00323     //
00324 
00325     // Storage for visited points
00326     labelHashSet visitedPoint(4*regionFaces.size());
00327 
00328     // Insert fence points so we don't visit them
00329     for
00330     (
00331         labelHashSet::const_iterator iter = fencePoints.begin();
00332         iter != fencePoints.end();
00333         ++iter
00334     )
00335     {
00336         visitedPoint.insert(iter.key());
00337     }
00338 
00339     labelHashSet visitedEdges(2*fencePoints.size());
00340 
00341 
00342     if (debug)
00343     {
00344         Info<< "Excluding visit of points:" << visitedPoint << endl;
00345     }
00346 
00347     for
00348     (
00349         labelHashSet::const_iterator iter = regionFaces.begin();
00350         iter != regionFaces.end();
00351         ++iter
00352     )
00353     {
00354         label faceI = iter.key();
00355 
00356         // Get side of face.
00357         label cellI;
00358 
00359         if (sideOwner_.found(faceI))
00360         {
00361             cellI = mesh.faceOwner()[faceI];
00362         }
00363         else
00364         {
00365             cellI = mesh.faceNeighbour()[faceI];
00366         }
00367 
00368         // Find starting point and edge on face.
00369         const labelList& fEdges = mesh.faceEdges()[faceI];
00370 
00371         forAll(fEdges, fEdgeI)
00372         {
00373             label edgeI = fEdges[fEdgeI];
00374 
00375             // Get the face 'perpendicular' to faceI on region.
00376             label otherFaceI = otherFace(mesh, cellI, faceI, edgeI);
00377 
00378             // Edge 
00379             const edge& e = mesh.edges()[edgeI];
00380 
00381             if (!visitedPoint.found(e.start()))
00382             {
00383                 Info<< "Determining visibility from point " << e.start()
00384                     << endl;
00385 
00386                 visitedPoint.insert(e.start());
00387 
00388                 //edgeI = otherEdge(mesh, otherFaceI, edgeI, e.start());
00389 
00390                 walkPointConnectedFaces
00391                 (
00392                     mesh,
00393                     regionEdges,
00394                     e.start(),
00395                     otherFaceI,
00396                     edgeI,
00397                     visitedEdges
00398                 );
00399             }
00400             if (!visitedPoint.found(e.end()))
00401             {
00402                 Info<< "Determining visibility from point " << e.end()
00403                     << endl;
00404 
00405                 visitedPoint.insert(e.end());
00406 
00407                 //edgeI = otherEdge(mesh, otherFaceI, edgeI, e.end());
00408 
00409                 walkPointConnectedFaces
00410                 (
00411                     mesh,
00412                     regionEdges,
00413                     e.end(),
00414                     otherFaceI,
00415                     edgeI,
00416                     visitedEdges
00417                 );
00418             }
00419         }
00420     }
00421 }
00422 
00423 
00424 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00425 
00426 // Construct from components
00427 Foam::regionSide::regionSide
00428 (
00429     const primitiveMesh& mesh,
00430     const labelHashSet& region,         // faces of region
00431     const labelHashSet& fenceEdges,     // outside edges
00432     const label startCellI,
00433     const label startFaceI
00434 )
00435 :
00436     sideOwner_(region.size()),
00437     insidePointFaces_(region.size())
00438 {
00439     // Storage for visited faces
00440     labelHashSet visitedFace(region.size());
00441 
00442     // Visit all faces on this side of region.
00443     // Mark for each face whether owner (or neighbour) has been visited.
00444     // Sets sideOwner_
00445     visitConnectedFaces
00446     (
00447         mesh,
00448         region,
00449         fenceEdges,
00450         startCellI,
00451         startFaceI,
00452         visitedFace
00453     );
00454 
00455     //
00456     // Visit all internal points on region and mark faces visitable from these.
00457     // Sets insidePointFaces_.
00458     //
00459 
00460     labelHashSet fencePoints(fenceEdges.size());
00461 
00462     for
00463     (
00464         labelHashSet::const_iterator iter = fenceEdges.begin();
00465         iter != fenceEdges.end();
00466         ++iter
00467     )
00468     {
00469         const edge& e = mesh.edges()[iter.key()];
00470 
00471         fencePoints.insert(e.start());
00472         fencePoints.insert(e.end());
00473     }
00474 
00475     walkAllPointConnectedFaces(mesh, region, fencePoints);
00476 }
00477 
00478 
00479 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines