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

surfaceToPatch.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 Application
00025     surfaceToPatch
00026 
00027 Description
00028     Reads surface and applies surface regioning to a mesh.
00029 
00030     Uses boundaryMesh to do the hard work.
00031 
00032 Usage
00033 
00034     - surfaceToPatch [OPTIONS] <Foam surface file>
00035 
00036     @param <Foam surface file> \n
00037     @todo Detailed description of argument.
00038 
00039     @param -tol <fraction of mesh size>\n
00040     Search toleracne.
00041 
00042     @param -faceSet <faceSet name>\n
00043     Only apply to named face set.
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 "boundaryMesh.H"
00062 #include <OpenFOAM/polyMesh.H>
00063 #include <meshTools/faceSet.H>
00064 #include <dynamicMesh/polyTopoChange.H>
00065 #include <dynamicMesh/polyModifyFace.H>
00066 #include <OpenFOAM/globalMeshData.H>
00067 
00068 using namespace Foam;
00069 
00070 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00071 
00072 // Adds empty patch if not yet there. Returns patchID.
00073 label addPatch(polyMesh& mesh, const word& patchName)
00074 {
00075     label patchI = mesh.boundaryMesh().findPatchID(patchName);
00076 
00077     if (patchI == -1)
00078     {
00079         const polyBoundaryMesh& patches = mesh.boundaryMesh();
00080 
00081         List<polyPatch*> newPatches(patches.size() + 1);
00082 
00083         patchI = 0;
00084 
00085         // Copy all old patches
00086         forAll(patches, i)
00087         {
00088             const polyPatch& pp = patches[i];
00089 
00090             newPatches[patchI] =
00091                 pp.clone
00092                 (
00093                     patches,
00094                     patchI,
00095                     pp.size(),
00096                     pp.start()
00097                 ).ptr();
00098 
00099             patchI++;
00100         }
00101 
00102         // Add zero-sized patch
00103         newPatches[patchI] =
00104             new polyPatch
00105             (
00106                 patchName,
00107                 0,
00108                 mesh.nFaces(),
00109                 patchI,
00110                 patches
00111             );
00112 
00113         mesh.removeBoundary();
00114         mesh.addPatches(newPatches);
00115 
00116         Pout<< "Created patch " << patchName << " at " << patchI << endl;
00117     }
00118     else
00119     {
00120         Pout<< "Reusing patch " << patchName << " at " << patchI << endl;
00121     }
00122 
00123     return patchI;
00124 }
00125 
00126 
00127 // Repatch single face. Return true if patch changed.
00128 bool repatchFace
00129 (
00130     const polyMesh& mesh,
00131     const boundaryMesh& bMesh,
00132     const labelList& nearest,
00133     const labelList& surfToMeshPatch,
00134     const label faceI,
00135     polyTopoChange& meshMod
00136 )
00137 {
00138     bool changed = false;
00139 
00140     label bFaceI = faceI - mesh.nInternalFaces();
00141 
00142     if (nearest[bFaceI] != -1)
00143     {
00144         // Use boundary mesh one.
00145         label bMeshPatchID = bMesh.whichPatch(nearest[bFaceI]);
00146 
00147         label patchID = surfToMeshPatch[bMeshPatchID];
00148 
00149         if (patchID != mesh.boundaryMesh().whichPatch(faceI))
00150         {
00151             label own = mesh.faceOwner()[faceI];
00152 
00153             label zoneID = mesh.faceZones().whichZone(faceI);
00154 
00155             bool zoneFlip = false;
00156 
00157             if (zoneID >= 0)
00158             {
00159                 const faceZone& fZone = mesh.faceZones()[zoneID];
00160 
00161                 zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
00162             }
00163 
00164             meshMod.setAction
00165             (
00166                 polyModifyFace
00167                 (
00168                     mesh.faces()[faceI],// modified face
00169                     faceI,              // label of face being modified
00170                     own,                // owner
00171                     -1,                 // neighbour
00172                     false,              // face flip
00173                     patchID,            // patch for face
00174                     false,              // remove from zone
00175                     zoneID,             // zone for face
00176                     zoneFlip            // face flip in zone
00177                 )
00178             );
00179 
00180             changed = true;
00181         }
00182     }
00183     else
00184     {
00185         changed = false;
00186     }
00187     return changed;
00188 }
00189 
00190 
00191 // Main program:
00192 
00193 int main(int argc, char *argv[])
00194 {
00195     argList::noParallel();
00196     argList::validArgs.append("surface file");
00197     argList::validOptions.insert("faceSet", "faceSet name");
00198     argList::validOptions.insert("tol", "fraction of mesh size");
00199 
00200 #   include <OpenFOAM/setRootCase.H>
00201 #   include <OpenFOAM/createTime.H>
00202 #   include <OpenFOAM/createPolyMesh.H>
00203 
00204     fileName surfName(args.additionalArgs()[0]);
00205 
00206     Info<< "Reading surface from " << surfName << " ..." << endl;
00207 
00208     bool readSet = args.optionFound("faceSet");
00209     word setName;
00210 
00211     if (readSet)
00212     {
00213         setName = args.option("faceSet");
00214 
00215         Info<< "Repatching only the faces in faceSet " << setName
00216             << " according to nearest surface triangle ..." << endl;
00217     }
00218     else
00219     {
00220         Info<< "Patching all boundary faces according to nearest surface"
00221             << " triangle ..." << endl;
00222     }
00223 
00224     scalar searchTol = 1E-3;
00225     args.optionReadIfPresent("tol", searchTol);
00226 
00227     // Get search box. Anything not within this box will not be considered.
00228     const boundBox& meshBb = mesh.globalData().bb();
00229     const vector searchSpan = searchTol * meshBb.span();
00230 
00231     Info<< "All boundary faces further away than " << searchTol
00232         << " of mesh bounding box " << meshBb
00233         << " will keep their patch label ..." << endl;
00234 
00235 
00236     Info<< "Before patching:" << nl
00237         << "    patch\tsize" << endl;
00238 
00239     forAll(mesh.boundaryMesh(), patchI)
00240     {
00241         Info<< "    " << mesh.boundaryMesh()[patchI].name() << '\t'
00242             << mesh.boundaryMesh()[patchI].size() << endl;
00243     }
00244     Info<< endl;
00245 
00246 
00247 
00248     boundaryMesh bMesh;
00249 
00250     // Load in the surface.
00251     bMesh.readTriSurface(surfName);
00252 
00253     // Add all the boundaryMesh patches to the mesh.
00254     const PtrList<boundaryPatch>& bPatches = bMesh.patches();
00255 
00256     // Map from surface patch ( = boundaryMesh patch) to polyMesh patch
00257     labelList patchMap(bPatches.size());
00258 
00259     forAll(bPatches, i)
00260     {
00261         patchMap[i] = addPatch(mesh, bPatches[i].name());
00262     }
00263 
00264     // Obtain nearest face in bMesh for each boundary face in mesh that
00265     // is within search span.
00266     // Note: should only determine for faceSet if working with that.
00267     labelList nearest(bMesh.getNearest(mesh, searchSpan));
00268 
00269     {
00270         // Dump unmatched faces to faceSet for debugging.
00271         faceSet unmatchedFaces(mesh, "unmatchedFaces", nearest.size()/100);
00272 
00273         forAll(nearest, bFaceI)
00274         {
00275             if (nearest[bFaceI] == -1)
00276             {
00277                 unmatchedFaces.insert(mesh.nInternalFaces() + bFaceI);
00278             }
00279         }
00280 
00281         Pout<< "Writing all " << unmatchedFaces.size()
00282             << " unmatched faces to faceSet "
00283             << unmatchedFaces.name()
00284             << endl;
00285 
00286         unmatchedFaces.write();
00287     }
00288 
00289 
00290     polyTopoChange meshMod(mesh);
00291 
00292     label nChanged = 0;
00293 
00294     if (readSet)
00295     {
00296         faceSet faceLabels(mesh, setName);
00297         Info<< "Read " << faceLabels.size() << " faces to repatch ..." << endl;
00298 
00299         forAllConstIter(faceSet, faceLabels, iter)
00300         {
00301             label faceI = iter.key();
00302 
00303             if (repatchFace(mesh, bMesh, nearest, patchMap, faceI, meshMod))
00304             {
00305                 nChanged++;
00306             }
00307         }
00308     }
00309     else
00310     {
00311         forAll(nearest, bFaceI)
00312         {
00313             label faceI = mesh.nInternalFaces() + bFaceI;
00314 
00315             if (repatchFace(mesh, bMesh, nearest, patchMap, faceI, meshMod))
00316             {
00317                 nChanged++;
00318             }
00319         }
00320     }
00321 
00322     Pout<< "Changed " << nChanged << " boundary faces." << nl << endl;
00323 
00324     if (nChanged > 0)
00325     {
00326         meshMod.changeMesh(mesh, false);
00327 
00328         Info<< "After patching:" << nl
00329             << "    patch\tsize" << endl;
00330 
00331         forAll(mesh.boundaryMesh(), patchI)
00332         {
00333             Info<< "    " << mesh.boundaryMesh()[patchI].name() << '\t'
00334                 << mesh.boundaryMesh()[patchI].size() << endl;
00335         }
00336         Info<< endl;
00337 
00338 
00339         runTime++;
00340 
00341         // Write resulting mesh
00342         Info << "Writing modified mesh to time " << runTime.value() << endl;
00343         mesh.write();
00344     }
00345 
00346 
00347     Info<< "End\n" << endl;
00348 
00349     return 0;
00350 }
00351 
00352 
00353 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines