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

combinePatchFaces.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     combinePatchFaces
00026 
00027 Description
00028     Checks for multiple patch faces on same cell and combines them.
00029 
00030     These result from e.g. refined neighbouring cells getting removed, leaving
00031     4 exposed faces with same owner.
00032 
00033     Rules for merging:
00034     - only boundary faces (since multiple internal faces between two cells
00035       not allowed anyway)
00036     - faces have to have same owner
00037     - faces have to be connected via edge which are not features (so angle
00038       between them < feature angle)
00039     - outside of faces has to be single loop
00040     - outside of face should not be (or just slightly) concave (so angle
00041       between consecutive edges < concaveangle
00042 
00043     E.g. to allow all faces on same patch to be merged:
00044 
00045         combinePatchFaces .. cavity 180 -concaveAngle 90
00046 
00047 Usage
00048 
00049     - combinePatchFaces [OPTIONS] <feature angle [0..180]>
00050 
00051     @param <feature angle [0..180]> \n
00052     @todo Detailed description of argument.
00053 
00054     @param -snapMesh \n
00055     Remove loose points on edges.
00056 
00057     @param -concaveAngle <angle [0..180]>\n
00058     Maximum allowed concave angle.
00059 
00060     @param -overwrite \n
00061     Overwrite existing data.
00062 
00063     @param -case <dir>\n
00064     Case directory.
00065 
00066     @param -parallel \n
00067     Run in parallel.
00068 
00069     @param -help \n
00070     Display help message.
00071 
00072     @param -doc \n
00073     Display Doxygen API documentation page for this application.
00074 
00075     @param -srcDoc \n
00076     Display Doxygen source documentation page for this application.
00077 
00078 \*---------------------------------------------------------------------------*/
00079 
00080 #include <OpenFOAM/PstreamReduceOps.H>
00081 #include <OpenFOAM/argList.H>
00082 #include <OpenFOAM/Time.H>
00083 #include <dynamicMesh/polyTopoChange.H>
00084 #include <dynamicMesh/polyModifyFace.H>
00085 #include <dynamicMesh/polyAddFace.H>
00086 #include <dynamicMesh/combineFaces.H>
00087 #include <dynamicMesh/removePoints.H>
00088 #include <OpenFOAM/polyMesh.H>
00089 #include <OpenFOAM/mapPolyMesh.H>
00090 #include <OpenFOAM/mathematicalConstants.H>
00091 
00092 using namespace Foam;
00093 
00094 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00095 
00096 // Sin of angle between two consecutive edges on a face. If sin(angle) larger
00097 // than this the face will be considered concave.
00098 const scalar defaultConcaveAngle = 30;
00099 
00100 
00101 // Same check as snapMesh
00102 void checkSnapMesh
00103 (
00104     const Time& runTime,
00105     const polyMesh& mesh,
00106     labelHashSet& wrongFaces
00107 )
00108 {
00109     IOdictionary snapDict
00110     (
00111         IOobject
00112         (
00113             "snapMeshDict",
00114             runTime.system(),
00115             mesh,
00116             IOobject::MUST_READ,
00117             IOobject::NO_WRITE
00118         )
00119     );
00120 
00121     // Max nonorthogonality allowed
00122     scalar maxNonOrtho(readScalar(snapDict.lookup("maxNonOrtho")));
00123     // Max concaveness allowed.
00124     scalar maxConcave(readScalar(snapDict.lookup("maxConcave")));
00125     // Min volume allowed (factor of minimum cellVolume)
00126     scalar relMinVol(readScalar(snapDict.lookup("minVol")));
00127     const scalar minCellVol = min(mesh.cellVolumes());
00128     const scalar minPyrVol = relMinVol*minCellVol;
00129     // Min area
00130     scalar minArea(readScalar(snapDict.lookup("minArea")));
00131 
00132     if (maxNonOrtho < 180.0-SMALL)
00133     {
00134         Pout<< "Checking non orthogonality" << endl;
00135 
00136         label nOldSize = wrongFaces.size();
00137         mesh.setNonOrthThreshold(maxNonOrtho);
00138         mesh.checkFaceOrthogonality(false, &wrongFaces);
00139 
00140         Pout<< "Detected " << wrongFaces.size() - nOldSize
00141             << " faces with non-orthogonality > " << maxNonOrtho << " degrees"
00142             << endl;
00143     }
00144 
00145     if (minPyrVol > -GREAT)
00146     {
00147         Pout<< "Checking face pyramids" << endl;
00148 
00149         label nOldSize = wrongFaces.size();
00150         mesh.checkFacePyramids(false, minPyrVol, &wrongFaces);
00151         Pout<< "Detected additional " << wrongFaces.size() - nOldSize
00152             << " faces with illegal face pyramids" << endl;
00153     }
00154 
00155     if (maxConcave < 180.0-SMALL)
00156     {
00157         Pout<< "Checking face angles" << endl;
00158 
00159         label nOldSize = wrongFaces.size();
00160         mesh.checkFaceAngles(false, maxConcave, &wrongFaces);
00161         Pout<< "Detected additional " << wrongFaces.size() - nOldSize
00162             << " faces with concavity > " << maxConcave << " degrees"
00163             << endl;
00164     }
00165 
00166     if (minArea > -SMALL)
00167     {
00168         Pout<< "Checking face areas" << endl;
00169 
00170         label nOldSize = wrongFaces.size();
00171 
00172         const scalarField magFaceAreas = mag(mesh.faceAreas());
00173 
00174         forAll(magFaceAreas, faceI)
00175         {
00176             if (magFaceAreas[faceI] < minArea)
00177             {
00178                 wrongFaces.insert(faceI);
00179             }
00180         }
00181         Pout<< "Detected additional " << wrongFaces.size() - nOldSize
00182             << " faces with area < " << minArea << " m^2" << endl;
00183     }
00184 }
00185 
00186 
00187 // Merge faces on the same patch (usually from exposing refinement)
00188 // Can undo merges if these cause problems.
00189 label mergePatchFaces
00190 (
00191     const scalar minCos,
00192     const scalar concaveSin,
00193     const bool snapMeshDict,
00194     const Time& runTime,
00195     polyMesh& mesh
00196 )
00197 {
00198     // Patch face merging engine
00199     combineFaces faceCombiner(mesh);
00200 
00201     // Get all sets of faces that can be merged
00202     labelListList allFaceSets(faceCombiner.getMergeSets(minCos, concaveSin));
00203 
00204     label nFaceSets = returnReduce(allFaceSets.size(), sumOp<label>());
00205 
00206     Info<< "Merging " << nFaceSets << " sets of faces." << endl;
00207 
00208     if (nFaceSets > 0)
00209     {
00210         // Store the faces of the face sets
00211         List<faceList> allFaceSetsFaces(allFaceSets.size());
00212         forAll(allFaceSets, setI)
00213         {
00214             allFaceSetsFaces[setI] = UIndirectList<face>
00215             (
00216                 mesh.faces(),
00217                 allFaceSets[setI]
00218             );
00219         }
00220 
00221         autoPtr<mapPolyMesh> map;
00222         {
00223             // Topology changes container
00224             polyTopoChange meshMod(mesh);
00225 
00226             // Merge all faces of a set into the first face of the set.
00227             faceCombiner.setRefinement(allFaceSets, meshMod);
00228 
00229             // Change the mesh (no inflation)
00230             map = meshMod.changeMesh(mesh, false, true);
00231 
00232             // Update fields
00233             mesh.updateMesh(map);
00234 
00235             // Move mesh (since morphing does not do this)
00236             if (map().hasMotionPoints())
00237             {
00238                 mesh.movePoints(map().preMotionPoints());
00239             }
00240             else
00241             {
00242                 // Delete mesh volumes. No other way to do this?
00243                 mesh.clearOut();
00244             }
00245         }
00246 
00247 
00248         // Check for errors and undo
00249         // ~~~~~~~~~~~~~~~~~~~~~~~~~
00250 
00251         // Faces in error.
00252         labelHashSet errorFaces;
00253 
00254         if (snapMeshDict)
00255         {
00256             checkSnapMesh(runTime, mesh, errorFaces);
00257         }
00258         else
00259         {
00260             mesh.checkFacePyramids(false, -SMALL, &errorFaces);
00261         }
00262 
00263         // Sets where the master is in error
00264         labelHashSet errorSets;
00265 
00266         forAll(allFaceSets, setI)
00267         {
00268             label newMasterI = map().reverseFaceMap()[allFaceSets[setI][0]];
00269 
00270             if (errorFaces.found(newMasterI))
00271             {
00272                 errorSets.insert(setI);
00273             }
00274         }
00275         label nErrorSets = returnReduce(errorSets.size(), sumOp<label>());
00276 
00277         Info<< "Detected " << nErrorSets
00278             << " error faces on boundaries that have been merged."
00279             << " These will be restored to their original faces."
00280             << endl;
00281 
00282         if (nErrorSets > 0)
00283         {
00284             // Renumber stored faces to new vertex numbering.
00285             forAllConstIter(labelHashSet, errorSets, iter)
00286             {
00287                 label setI = iter.key();
00288 
00289                 faceList& setFaceVerts = allFaceSetsFaces[setI];
00290 
00291                 forAll(setFaceVerts, i)
00292                 {
00293                     inplaceRenumber(map().reversePointMap(), setFaceVerts[i]);
00294 
00295                     // Debug: check that all points are still there.
00296                     forAll(setFaceVerts[i], j)
00297                     {
00298                         label newVertI = setFaceVerts[i][j];
00299 
00300                         if (newVertI < 0)
00301                         {
00302                             FatalErrorIn("mergePatchFaces")
00303                                 << "In set:" << setI << " old face labels:"
00304                                 << allFaceSets[setI] << " new face vertices:"
00305                                 << setFaceVerts[i] << " are unmapped vertices!"
00306                                 << abort(FatalError);
00307                         }
00308                     }
00309                 }
00310             }
00311 
00312 
00313             // Topology changes container
00314             polyTopoChange meshMod(mesh);
00315 
00316 
00317             // Restore faces
00318             forAllConstIter(labelHashSet, errorSets, iter)
00319             {
00320                 label setI = iter.key();
00321 
00322                 const labelList& setFaces = allFaceSets[setI];
00323                 const faceList& setFaceVerts = allFaceSetsFaces[setI];
00324 
00325                 label newMasterI = map().reverseFaceMap()[setFaces[0]];
00326 
00327                 // Restore. Get face properties.
00328 
00329                 label own = mesh.faceOwner()[newMasterI];
00330                 label zoneID = mesh.faceZones().whichZone(newMasterI);
00331                 bool zoneFlip = false;
00332                 if (zoneID >= 0)
00333                 {
00334                     const faceZone& fZone = mesh.faceZones()[zoneID];
00335                     zoneFlip = fZone.flipMap()[fZone.whichFace(newMasterI)];
00336                 }
00337                 label patchI = mesh.boundaryMesh().whichPatch(newMasterI);
00338 
00339 
00340                 Pout<< "Restoring new master face " << newMasterI
00341                     << " to vertices " << setFaceVerts[0] << endl;
00342 
00343                 // Modify the master face.
00344                 meshMod.setAction
00345                 (
00346                     polyModifyFace
00347                     (
00348                         setFaceVerts[0],                // original face
00349                         newMasterI,                     // label of face
00350                         own,                            // owner
00351                         -1,                             // neighbour
00352                         false,                          // face flip
00353                         patchI,                         // patch for face
00354                         false,                          // remove from zone
00355                         zoneID,                         // zone for face
00356                         zoneFlip                        // face flip in zone
00357                     )
00358                 );
00359 
00360 
00361                 // Add the previously removed faces
00362                 for (label i = 1; i < setFaces.size(); i++)
00363                 {
00364                     Pout<< "Restoring removed face " << setFaces[i]
00365                         << " with vertices " << setFaceVerts[i] << endl;
00366 
00367                     meshMod.setAction
00368                     (
00369                         polyAddFace
00370                         (
00371                             setFaceVerts[i],        // vertices
00372                             own,                    // owner,
00373                             -1,                     // neighbour,
00374                             -1,                     // masterPointID,
00375                             -1,                     // masterEdgeID,
00376                             newMasterI,             // masterFaceID,
00377                             false,                  // flipFaceFlux,
00378                             patchI,                 // patchID,
00379                             zoneID,                 // zoneID,
00380                             zoneFlip                // zoneFlip
00381                         )
00382                     );
00383                 }
00384             }
00385 
00386             // Change the mesh (no inflation)
00387             map = meshMod.changeMesh(mesh, false, true);
00388 
00389             // Update fields
00390             mesh.updateMesh(map);
00391 
00392             // Move mesh (since morphing does not do this)
00393             if (map().hasMotionPoints())
00394             {
00395                 mesh.movePoints(map().preMotionPoints());
00396             }
00397             else
00398             {
00399                 // Delete mesh volumes. No other way to do this?
00400                 mesh.clearOut();
00401             }
00402         }
00403     }
00404     else
00405     {
00406         Info<< "No faces merged ..." << endl;
00407     }
00408 
00409     return nFaceSets;
00410 }
00411 
00412 
00413 // Remove points not used by any face or points used by only two faces where
00414 // the edges are in line
00415 label mergeEdges(const scalar minCos, polyMesh& mesh)
00416 {
00417     Info<< "Merging all points on surface that" << nl
00418         << "- are used by only two boundary faces and" << nl
00419         << "- make an angle with a cosine of more than " << minCos
00420         << "." << nl << endl;
00421 
00422     // Point removal analysis engine
00423     removePoints pointRemover(mesh);
00424 
00425     // Count usage of points
00426     boolList pointCanBeDeleted;
00427     label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
00428 
00429     if (nRemove > 0)
00430     {
00431         Info<< "Removing " << nRemove
00432             << " straight edge points ..." << endl;
00433 
00434         // Topology changes container
00435         polyTopoChange meshMod(mesh);
00436 
00437         pointRemover.setRefinement(pointCanBeDeleted, meshMod);
00438 
00439         // Change the mesh (no inflation)
00440         autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
00441 
00442         // Update fields
00443         mesh.updateMesh(map);
00444 
00445         // Move mesh (since morphing does not do this)
00446         if (map().hasMotionPoints())
00447         {
00448             mesh.movePoints(map().preMotionPoints());
00449         }
00450         else
00451         {
00452             // Delete mesh volumes. No other way to do this?
00453             mesh.clearOut();
00454         }
00455     }
00456     else
00457     {
00458         Info<< "No straight edges simplified and no points removed ..." << endl;
00459     }
00460 
00461     return nRemove;
00462 }
00463 
00464 
00465 // Main program:
00466 
00467 int main(int argc, char *argv[])
00468 {
00469     argList::validArgs.append("feature angle [0..180]");
00470     argList::validOptions.insert("concaveAngle", "[0..180]");
00471     argList::validOptions.insert("snapMesh", "");
00472     argList::validOptions.insert("overwrite", "");
00473 
00474 #   include <OpenFOAM/setRootCase.H>
00475 #   include <OpenFOAM/createTime.H>
00476     runTime.functionObjects().off();
00477 #   include <OpenFOAM/createPolyMesh.H>
00478     const word oldInstance = mesh.pointsInstance();
00479 
00480     scalar featureAngle(readScalar(IStringStream(args.additionalArgs()[0])()));
00481 
00482     scalar minCos = Foam::cos(featureAngle*mathematicalConstant::pi/180.0);
00483 
00484     scalar concaveAngle = defaultConcaveAngle;
00485     args.optionReadIfPresent("concaveAngle", concaveAngle);
00486 
00487     scalar concaveSin = Foam::sin(concaveAngle*mathematicalConstant::pi/180.0);
00488 
00489     bool snapMeshDict = args.optionFound("snapMesh");
00490     bool overwrite = args.optionFound("overwrite");
00491 
00492     Info<< "Merging all faces of a cell" << nl
00493         << "    - which are on the same patch" << nl
00494         << "    - which make an angle < " << featureAngle << " degrees"
00495         << nl
00496         << "      (cos:" << minCos << ')' << nl
00497         << "    - even when resulting face becomes concave by more than "
00498         << concaveAngle << " degrees" << nl
00499         << "      (sin:" << concaveSin << ')' << nl
00500         << endl;
00501 
00502     if (!overwrite)
00503     {
00504         runTime++;
00505     }
00506 
00507     // Merge faces on same patch
00508     label nChanged = mergePatchFaces
00509     (
00510         minCos,
00511         concaveSin,
00512         snapMeshDict,
00513         runTime,
00514         mesh
00515     );
00516 
00517     // Merge points on straight edges and remove unused points
00518     if (snapMeshDict)
00519     {
00520         Info<< "Merging all 'loose' points on surface edges"
00521             << ", regardless of the angle they make." << endl;
00522 
00523         // Surface bnound to be used to extrude. Merge all loose points.
00524         nChanged += mergeEdges(-1, mesh);
00525     }
00526     else
00527     {
00528         nChanged += mergeEdges(minCos, mesh);
00529     }
00530 
00531     if (nChanged > 0)
00532     {
00533         if (overwrite)
00534         {
00535             mesh.setInstance(oldInstance);
00536         }
00537 
00538         Info<< "Writing morphed mesh to time " << runTime.timeName() << endl;
00539 
00540         mesh.write();
00541     }
00542     else
00543     {
00544         Info<< "Mesh unchanged." << endl;
00545     }
00546 
00547     Info << "End\n" << endl;
00548 
00549     return 0;
00550 }
00551 
00552 
00553 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines