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

selectCells.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-2011 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 Application
00025     selectCells
00026 
00027 Description
00028     Select cells in relation to surface.
00029 
00030     Divides cells into three sets:
00031     - cutCells : cells cut by surface or close to surface.
00032     - outside  : cells not in cutCells and reachable from set of
00033       user-defined points (outsidePoints)
00034     - inside   : same but not reachable.
00035 
00036     Finally the wanted sets are combined into a cellSet 'selected'. Apart
00037     from straightforward adding the contents there are a few extra rules to
00038     make sure that the surface of the 'outside' of the mesh is singly
00039     connected.
00040 
00041 Usage
00042 
00043     - selectCells [OPTIONS]
00044 
00045     @param -case <dir>\n
00046     Case directory.
00047 
00048     @param -help \n
00049     Display help message.
00050 
00051     @param -doc \n
00052     Display Doxygen API documentation page for this application.
00053 
00054     @param -srcDoc \n
00055     Display Doxygen source documentation page for this application.
00056 
00057 \*---------------------------------------------------------------------------*/
00058 
00059 #include <OpenFOAM/argList.H>
00060 #include <OpenFOAM/Time.H>
00061 #include <OpenFOAM/polyMesh.H>
00062 #include <OpenFOAM/IOdictionary.H>
00063 #include <meshTools/twoDPointCorrector.H>
00064 #include <OpenFOAM/OFstream.H>
00065 #include <meshTools/meshTools.H>
00066 
00067 #include <triSurface/triSurface.H>
00068 #include <meshTools/triSurfaceSearch.H>
00069 #include <meshTools/meshSearch.H>
00070 #include <meshTools/cellClassification.H>
00071 #include <meshTools/cellSet.H>
00072 #include <meshTools/cellInfo.H>
00073 #include <OpenFOAM/MeshWave.H>
00074 #include "edgeStats.H"
00075 #include <meshTools/treeDataTriSurface.H>
00076 #include <meshTools/indexedOctree.H>
00077 
00078 using namespace Foam;
00079 
00080 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00081 
00082 // cellType for cells included/not included in mesh.
00083 static const label MESH = cellClassification::INSIDE;
00084 static const label NONMESH = cellClassification::OUTSIDE;
00085 
00086 
00087 void writeSet(const cellSet& cells, const string& msg)
00088 {
00089     Info<< "Writing " << msg << " (" << cells.size() << ") to cellSet "
00090         << cells.instance()/cells.local()/cells.name()
00091         << endl << endl;
00092     cells.write();
00093 }
00094 
00095 
00096 void getType(const labelList& elems, const label type, labelHashSet& set)
00097 {
00098     forAll(elems, i)
00099     {
00100         if (elems[i] == type)
00101         {
00102             set.insert(i);
00103         }
00104     }
00105 }
00106 
00107 
00108 void cutBySurface
00109 (
00110     const polyMesh& mesh,
00111     const meshSearch& queryMesh,
00112     const triSurfaceSearch& querySurf,
00113 
00114     const pointField& outsidePts,
00115     const bool selectCut,
00116     const bool selectInside,
00117     const bool selectOutside,
00118     const scalar nearDist,
00119 
00120     cellClassification& cellType
00121 )
00122 {
00123     // Cut with surface and classify as inside/outside/cut
00124     cellType =
00125         cellClassification
00126         (
00127             mesh,
00128             queryMesh,
00129             querySurf,
00130             outsidePts
00131         );
00132 
00133     // Get inside/outside/cutCells cellSets.
00134     cellSet inside(mesh, "inside", mesh.nCells()/10);
00135     getType(cellType, cellClassification::INSIDE, inside);
00136     writeSet(inside, "inside cells");
00137 
00138     cellSet outside(mesh, "outside", mesh.nCells()/10);
00139     getType(cellType, cellClassification::OUTSIDE, outside);
00140     writeSet(outside, "outside cells");
00141 
00142     cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
00143     getType(cellType, cellClassification::CUT, cutCells);
00144     writeSet(cutCells, "cells cut by surface");
00145 
00146 
00147     // Change cellType to reflect selected part of mesh. Use
00148     // MESH to denote selected part, NONMESH for all
00149     // other cells.
00150     // Is a bit of a hack but allows us to reuse all the functionality
00151     // in cellClassification.
00152 
00153     forAll(cellType, cellI)
00154     {
00155         label cType = cellType[cellI];
00156 
00157         if (cType == cellClassification::CUT)
00158         {
00159             if (selectCut)
00160             {
00161                 cellType[cellI] = MESH;
00162             }
00163             else
00164             {
00165                 cellType[cellI] = NONMESH;
00166             }
00167         }
00168         else if (cType == cellClassification::INSIDE)
00169         {
00170             if (selectInside)
00171             {
00172                 cellType[cellI] = MESH;
00173             }
00174             else
00175             {
00176                 cellType[cellI] = NONMESH;
00177             }
00178         }
00179         else if (cType == cellClassification::OUTSIDE)
00180         {
00181             if (selectOutside)
00182             {
00183                 cellType[cellI] = MESH;
00184             }
00185             else
00186             {
00187                 cellType[cellI] = NONMESH;
00188             }
00189         }
00190         else
00191         {
00192             FatalErrorIn("cutBySurface")
00193                 << "Multiple mesh regions in original mesh" << endl
00194                 << "Please use splitMeshRegions to separate these"
00195                 << exit(FatalError);
00196         }
00197     }
00198 
00199 
00200     if (nearDist > 0)
00201     {
00202         Info<< "Removing cells with points closer than " << nearDist
00203             << " to the surface ..." << nl << endl;
00204 
00205         const pointField& pts = mesh.points();
00206         const indexedOctree<treeDataTriSurface>& tree = querySurf.tree();
00207 
00208         label nRemoved = 0;
00209 
00210         forAll(pts, pointI)
00211         {
00212             const point& pt = pts[pointI];
00213 
00214             pointIndexHit hitInfo = tree.findNearest(pt, sqr(nearDist));
00215 
00216             if (hitInfo.hit())
00217             {
00218                 const labelList& pCells = mesh.pointCells()[pointI];
00219 
00220                 forAll(pCells, i)
00221                 {
00222                     if (cellType[pCells[i]] != NONMESH)
00223                     {
00224                         cellType[pCells[i]] = NONMESH;
00225                         nRemoved++;
00226                     }
00227                 }
00228             }
00229         }
00230 
00231 //        tmp<pointField> tnearest = querySurf.calcNearest(pts);
00232 //        const pointField& nearest = tnearest();
00233 //
00234 //        label nRemoved = 0;
00235 //
00236 //        forAll(nearest, pointI)
00237 //        {
00238 //            if (mag(nearest[pointI] - pts[pointI]) < nearDist)
00239 //            {
00240 //                const labelList& pCells = mesh.pointCells()[pointI];
00241 //
00242 //                forAll(pCells, i)
00243 //                {
00244 //                    if (cellType[pCells[i]] != NONMESH)
00245 //                    {
00246 //                        cellType[pCells[i]] = NONMESH;
00247 //                        nRemoved++;
00248 //                    }
00249 //                }
00250 //            }
00251 //        }
00252 
00253         Info<< "Removed " << nRemoved << " cells since too close to surface"
00254             << nl << endl;
00255     }
00256 }
00257 
00258 
00259 
00260 // We're meshing the outside. Subset the currently selected mesh cells with the
00261 // ones reachable from the outsidepoints.
00262 label selectOutsideCells
00263 (
00264     const polyMesh& mesh,
00265     const meshSearch& queryMesh,
00266     const pointField& outsidePts,
00267     cellClassification& cellType
00268 )
00269 {
00270     //
00271     // Check all outsidePts and for all of them inside a mesh cell
00272     // collect the faces to start walking from
00273     //
00274 
00275     // Outside faces
00276     labelHashSet outsideFacesMap(outsidePts.size() * 6 * 2);
00277     DynamicList<label> outsideFaces(outsideFacesMap.size());
00278     // CellInfo on outside faces
00279     DynamicList<cellInfo> outsideFacesInfo(outsideFacesMap.size());
00280 
00281     // cellInfo for mesh cell
00282     const cellInfo meshInfo(MESH);
00283 
00284     forAll(outsidePts, outsidePtI)
00285     {
00286         // Find cell containing point. Linear search.
00287         label cellI = queryMesh.findCell(outsidePts[outsidePtI], -1, false);
00288 
00289         if (cellI != -1 && cellType[cellI] == MESH)
00290         {
00291             Info<< "Marking cell " << cellI << " containing outside point "
00292                 << outsidePts[outsidePtI] << " with type " << cellType[cellI]
00293                 << " ..." << endl;
00294 
00295             //
00296             // Mark this cell and its faces to start walking from
00297             //
00298 
00299             // Mark faces of cellI
00300             const labelList& cFaces = mesh.cells()[cellI];
00301             forAll(cFaces, i)
00302             {
00303                 label faceI = cFaces[i];
00304 
00305                 if (outsideFacesMap.insert(faceI))
00306                 {
00307                     outsideFaces.append(faceI);
00308                     outsideFacesInfo.append(meshInfo);
00309                 }
00310             }
00311         }
00312     }
00313 
00314     // Floodfill starting from outsideFaces (of type meshInfo)
00315     MeshWave<cellInfo> regionCalc
00316     (
00317         mesh,
00318         outsideFaces.shrink(),
00319         outsideFacesInfo.shrink(),
00320         mesh.nCells()  // max iterations
00321     );
00322 
00323     // Now regionCalc should hold info on cells that are reachable from
00324     // changedFaces. Use these to subset cellType
00325     const List<cellInfo>& allCellInfo = regionCalc.allCellInfo();
00326 
00327     label nChanged = 0;
00328 
00329     forAll(allCellInfo, cellI)
00330     {
00331         if (cellType[cellI] == MESH)
00332         {
00333             // Original cell was selected for meshing. Check if cell was
00334             // reached from outsidePoints
00335             if (allCellInfo[cellI].type() != MESH)
00336             {
00337                 cellType[cellI] = NONMESH;
00338                 nChanged++;
00339             }
00340         }
00341     }
00342 
00343     return nChanged;
00344 }
00345 
00346 
00347 // Main program:
00348 
00349 int main(int argc, char *argv[])
00350 {
00351     argList::noParallel();
00352 
00353 #   include <OpenFOAM/setRootCase.H>
00354 #   include <OpenFOAM/createTime.H>
00355 #   include <OpenFOAM/createPolyMesh.H>
00356 
00357     // Mesh edge statistics calculator
00358     edgeStats edgeCalc(mesh);
00359 
00360 
00361     IOdictionary refineDict
00362     (
00363         IOobject
00364         (
00365             "selectCellsDict",
00366             runTime.system(),
00367             mesh,
00368             IOobject::MUST_READ,
00369             IOobject::NO_WRITE
00370         )
00371     );
00372 
00373     fileName surfName(refineDict.lookup("surface"));
00374     pointField outsidePts(refineDict.lookup("outsidePoints"));
00375     bool useSurface(readBool(refineDict.lookup("useSurface")));
00376     bool selectCut(readBool(refineDict.lookup("selectCut")));
00377     bool selectInside(readBool(refineDict.lookup("selectInside")));
00378     bool selectOutside(readBool(refineDict.lookup("selectOutside")));
00379     scalar nearDist(readScalar(refineDict.lookup("nearDistance")));
00380 
00381 
00382     if (useSurface)
00383     {
00384         Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
00385             << "    cells cut by surface            : " << selectCut << nl
00386             << "    cells inside of surface         : " << selectInside << nl
00387             << "    cells outside of surface        : " << selectOutside << nl
00388             << "    cells with points further than  : " << nearDist << nl
00389             << endl;
00390     }
00391     else
00392     {
00393         Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
00394             << "    cells reachable from outsidePoints:" << selectOutside << nl
00395             << endl;
00396     }
00397 
00398     // Print edge stats on original mesh.
00399     (void)edgeCalc.minLen(Info);
00400 
00401     // Search engine on mesh. Face decomposition since faces might be warped.
00402     meshSearch queryMesh(mesh, true);
00403 
00404     // Check all 'outside' points
00405     forAll(outsidePts, outsideI)
00406     {
00407         const point& outsidePoint = outsidePts[outsideI];
00408 
00409         label cellI = queryMesh.findCell(outsidePoint, -1, false);
00410         if (returnReduce(cellI, maxOp<label>()) == -1)
00411         {
00412             FatalErrorIn(args.executable())
00413                 << "outsidePoint " << outsidePoint
00414                 << " is not inside any cell"
00415                 << exit(FatalError);
00416         }
00417     }
00418 
00419     // Cell status (compared to surface if provided): inside/outside/cut.
00420     // Start off from everything selected and cut later.
00421     cellClassification cellType
00422     (
00423         mesh,
00424         labelList
00425         (
00426             mesh.nCells(),
00427             cellClassification::MESH
00428         )
00429     );
00430 
00431 
00432     // Surface
00433     autoPtr<triSurface> surf(NULL);
00434     // Search engine on surface.
00435     autoPtr<triSurfaceSearch> querySurf(NULL);
00436 
00437     if (useSurface)
00438     {
00439         surf.reset(new triSurface(surfName));
00440 
00441         // Dump some stats
00442         surf().writeStats(Info);
00443 
00444         // Search engine on surface.
00445         querySurf.reset(new triSurfaceSearch(surf));
00446 
00447         // Set cellType[cellI] according to relation to surface
00448         cutBySurface
00449         (
00450             mesh,
00451             queryMesh,
00452             querySurf,
00453 
00454             outsidePts,
00455             selectCut,
00456             selectInside,
00457             selectOutside,
00458             nearDist,
00459 
00460             cellType
00461         );
00462     }
00463 
00464 
00465     // Now 'trim' all the corners from the mesh so meshing/surface extraction
00466     // becomes easier.
00467 
00468     label nHanging, nRegionEdges, nRegionPoints, nOutside;
00469 
00470     do
00471     {
00472         Info<< "Removing cells which after subsetting would have all points"
00473             << " on outside ..." << nl << endl;
00474 
00475         nHanging = cellType.fillHangingCells
00476         (
00477             MESH,       // meshType
00478             NONMESH,    // fill type
00479             mesh.nCells()
00480         );
00481 
00482 
00483         Info<< "Removing edges connecting cells unconnected by faces ..."
00484             << nl << endl;
00485 
00486         nRegionEdges = cellType.fillRegionEdges
00487         (
00488             MESH,       // meshType
00489             NONMESH,    // fill type
00490             mesh.nCells()
00491         ); 
00492 
00493 
00494         Info<< "Removing points connecting cells unconnected by faces ..."
00495             << nl << endl;
00496 
00497         nRegionPoints = cellType.fillRegionPoints
00498         (
00499             MESH,       // meshType
00500             NONMESH,    // fill type
00501             mesh.nCells()
00502         );
00503 
00504         nOutside = 0;
00505         if (selectOutside)
00506         {
00507             // Since we're selecting the cells reachable from outsidePoints
00508             // and the set might have changed, redo the outsideCells
00509             // calculation
00510             nOutside = selectOutsideCells
00511             (
00512                 mesh,
00513                 queryMesh,
00514                 outsidePts,
00515                 cellType
00516             );
00517         }
00518     } while
00519     (
00520         nHanging != 0
00521      || nRegionEdges != 0
00522      || nRegionPoints != 0
00523      || nOutside != 0
00524     );
00525 
00526     cellSet selectedCells(mesh, "selected", mesh.nCells()/10);
00527     getType(cellType, MESH, selectedCells);
00528 
00529     writeSet(selectedCells, "cells selected for meshing");
00530 
00531 
00532     Info << "End\n" << endl;
00533 
00534     return 0;
00535 }
00536 
00537 
00538 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines