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

modifyMesh.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     modifyMesh
00026 
00027 Description
00028     Manipulates mesh elements.
00029 
00030     Actions are:
00031         (boundary)points:
00032             - move
00033 
00034         (boundary)edges:
00035             - split and move introduced point
00036 
00037         (boundary)faces:
00038             - split(triangulate) and move introduced point
00039 
00040         edges:
00041             - collapse
00042 
00043         cells:
00044             - split into polygonal base pyramids around newly introduced mid
00045               point
00046 
00047     Is a bit of a loose collection of mesh change drivers.
00048 
00049 Usage
00050 
00051     - modifyMesh [OPTIONS]
00052 
00053     @param -overwrite \n
00054     Overwrite existing data.
00055 
00056     @param -case <dir>\n
00057     Case directory.
00058 
00059     @param -parallel \n
00060     Run in parallel.
00061 
00062     @param -help \n
00063     Display help message.
00064 
00065     @param -doc \n
00066     Display Doxygen API documentation page for this application.
00067 
00068     @param -srcDoc \n
00069     Display Doxygen source documentation page for this application.
00070 
00071 \*---------------------------------------------------------------------------*/
00072 
00073 #include <OpenFOAM/argList.H>
00074 #include <OpenFOAM/Time.H>
00075 #include <OpenFOAM/polyMesh.H>
00076 #include <dynamicMesh/polyTopoChange.H>
00077 #include <OpenFOAM/mapPolyMesh.H>
00078 #include <dynamicMesh/boundaryCutter.H>
00079 #include "cellSplitter.H"
00080 #include <dynamicMesh/edgeCollapser.H>
00081 #include <meshTools/meshTools.H>
00082 #include <OpenFOAM/Pair.H>
00083 
00084 using namespace Foam;
00085 
00086 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00087 
00088 // Locate point on patch. Returns (mesh) point label.
00089 label findPoint(const primitivePatch& pp, const point& nearPoint)
00090 {
00091     const pointField& points = pp.points();
00092     const labelList& meshPoints = pp.meshPoints();
00093 
00094     // Find nearest and next nearest
00095     scalar minDistSqr = GREAT;
00096     label minI = -1;
00097 
00098     scalar almostMinDistSqr = GREAT;
00099     label almostMinI = -1;
00100 
00101     forAll(meshPoints, i)
00102     {
00103         label pointI = meshPoints[i];
00104 
00105         scalar distSqr = magSqr(nearPoint - points[pointI]);
00106 
00107         if (distSqr < minDistSqr)
00108         {
00109             almostMinDistSqr = minDistSqr;
00110             almostMinI = minI;
00111 
00112             minDistSqr = distSqr;
00113             minI = pointI;
00114         }
00115         else if (distSqr < almostMinDistSqr)
00116         {
00117             almostMinDistSqr = distSqr;
00118             almostMinI = pointI;
00119         }
00120     }
00121 
00122 
00123     // Decide if nearPoint unique enough.
00124     Info<< "Found to point " << nearPoint << nl
00125         << "    nearest point      : " << minI
00126         << " distance " <<  Foam::sqrt(minDistSqr)
00127         << " at " << points[minI] << nl
00128         << "    next nearest point : " << almostMinI
00129         << " distance " <<  Foam::sqrt(almostMinDistSqr)
00130         << " at " << points[almostMinI] << endl;
00131 
00132     if (almostMinDistSqr < 4*minDistSqr)
00133     {
00134         Info<< "Next nearest too close to nearest. Aborting" << endl;
00135 
00136         return -1;
00137     }
00138     else
00139     {
00140         return minI;
00141     }
00142 }
00143 
00144 
00145 // Locate edge on patch. Return mesh edge label.
00146 label findEdge
00147 (
00148     const primitiveMesh& mesh,
00149     const primitivePatch& pp,
00150     const point& nearPoint
00151 )
00152 {
00153     const pointField& localPoints = pp.localPoints();
00154     const pointField& points = pp.points();
00155     const labelList& meshPoints = pp.meshPoints();
00156     const edgeList& edges = pp.edges();
00157 
00158     // Find nearest and next nearest
00159     scalar minDist = GREAT;
00160     label minI = -1;
00161 
00162     scalar almostMinDist = GREAT;
00163     label almostMinI = -1;
00164 
00165     forAll(edges, edgeI)
00166     {
00167         const edge& e = edges[edgeI];
00168 
00169         pointHit pHit(e.line(localPoints).nearestDist(nearPoint));
00170 
00171         if (pHit.hit())
00172         {
00173             if (pHit.distance() < minDist)
00174             {
00175                 almostMinDist = minDist;
00176                 almostMinI = minI;
00177 
00178                 minDist = pHit.distance();
00179                 minI = meshTools::findEdge
00180                 (
00181                     mesh,
00182                     meshPoints[e[0]],
00183                     meshPoints[e[1]]
00184                 );
00185             }
00186             else if (pHit.distance() < almostMinDist)
00187             {
00188                 almostMinDist = pHit.distance();
00189                 almostMinI = meshTools::findEdge
00190                 (
00191                     mesh,
00192                     meshPoints[e[0]],
00193                     meshPoints[e[1]]
00194                 );
00195             }
00196         }
00197     }
00198 
00199     if (minI == -1)
00200     {
00201         Info<< "Did not find edge close to point " << nearPoint << endl;
00202 
00203         return -1;
00204     }
00205 
00206 
00207     // Decide if nearPoint unique enough.
00208     Info<< "Found to point " << nearPoint << nl
00209         << "    nearest edge      : " << minI
00210         << " distance " << minDist << " endpoints "
00211         << mesh.edges()[minI].line(points) << nl
00212         << "    next nearest edge : " << almostMinI
00213         << " distance " << almostMinDist << " endpoints "
00214         << mesh.edges()[almostMinI].line(points) << nl
00215         << endl;
00216 
00217     if (almostMinDist < 2*minDist)
00218     {
00219         Info<< "Next nearest too close to nearest. Aborting" << endl;
00220 
00221         return -1;
00222     }
00223     else
00224     {
00225         return minI;
00226     }
00227 }
00228 
00229 
00230 // Find face on patch. Return mesh face label.
00231 label findFace
00232 (
00233     const primitiveMesh& mesh,
00234     const primitivePatch& pp,
00235     const point& nearPoint
00236 )
00237 {
00238     const pointField& points = pp.points();
00239 
00240     // Find nearest and next nearest
00241     scalar minDist = GREAT;
00242     label minI = -1;
00243 
00244     scalar almostMinDist = GREAT;
00245     label almostMinI = -1;
00246 
00247     forAll(pp, patchFaceI)
00248     {
00249         pointHit pHit(pp[patchFaceI].nearestPoint(nearPoint, points));
00250 
00251         if (pHit.hit())
00252         {
00253             if (pHit.distance() < minDist)
00254             {
00255                 almostMinDist = minDist;
00256                 almostMinI = minI;
00257 
00258                 minDist = pHit.distance();
00259                 minI = patchFaceI + mesh.nInternalFaces();
00260             }
00261             else if (pHit.distance() < almostMinDist)
00262             {
00263                 almostMinDist = pHit.distance();
00264                 almostMinI = patchFaceI + mesh.nInternalFaces();
00265             }
00266         }
00267     }
00268 
00269     if (minI == -1)
00270     {
00271         Info<< "Did not find face close to point " << nearPoint << endl;
00272 
00273         return -1;
00274     }
00275 
00276 
00277     // Decide if nearPoint unique enough.
00278     Info<< "Found to point " << nearPoint << nl
00279         << "    nearest face      : " << minI
00280         << " distance " << minDist
00281         << " to face centre " << mesh.faceCentres()[minI] << nl
00282         << "    next nearest face : " << almostMinI
00283         << " distance " << almostMinDist
00284         << " to face centre " << mesh.faceCentres()[almostMinI] << nl
00285         << endl;
00286 
00287     if (almostMinDist < 2*minDist)
00288     {
00289         Info<< "Next nearest too close to nearest. Aborting" << endl;
00290 
00291         return -1;
00292     }
00293     else
00294     {
00295         return minI;
00296     }
00297 }
00298 
00299 
00300 // Find cell with cell centre close to given point.
00301 label findCell(const primitiveMesh& mesh, const point& nearPoint)
00302 {
00303     label cellI = mesh.findCell(nearPoint);
00304 
00305     if (cellI != -1)
00306     {
00307         scalar distToCcSqr = magSqr(nearPoint - mesh.cellCentres()[cellI]);
00308 
00309         const labelList& cPoints = mesh.cellPoints()[cellI];
00310 
00311         label minI = -1;
00312         scalar minDistSqr = GREAT;
00313 
00314         forAll(cPoints, i)
00315         {
00316             label pointI = cPoints[i];
00317 
00318             scalar distSqr = magSqr(nearPoint - mesh.points()[pointI]);
00319 
00320             if (distSqr < minDistSqr)
00321             {
00322                 minDistSqr = distSqr;
00323                 minI = pointI;
00324             }
00325         }
00326 
00327         // Decide if nearPoint unique enough.
00328         Info<< "Found to point " << nearPoint << nl
00329             << "    nearest cell       : " << cellI
00330             << " distance " << Foam::sqrt(distToCcSqr)
00331             << " to cell centre " << mesh.cellCentres()[cellI] << nl
00332             << "    nearest mesh point : " << minI
00333             << " distance " << Foam::sqrt(minDistSqr)
00334             << " to " << mesh.points()[minI] << nl
00335             << endl;
00336 
00337         if (minDistSqr < 4*distToCcSqr)
00338         {
00339             Info<< "Mesh point too close to nearest cell centre. Aborting"
00340                 << endl;
00341 
00342             cellI = -1;
00343         }
00344     }
00345 
00346     return cellI;
00347 }
00348 
00349 
00350 
00351 // Main program:
00352 
00353 int main(int argc, char *argv[])
00354 {
00355     argList::validOptions.insert("overwrite", "");
00356 
00357 #   include <OpenFOAM/setRootCase.H>
00358 #   include <OpenFOAM/createTime.H>
00359     runTime.functionObjects().off();
00360 #   include <OpenFOAM/createPolyMesh.H>
00361     const word oldInstance = mesh.pointsInstance();
00362 
00363     bool overwrite = args.optionFound("overwrite");
00364 
00365     Info<< "Reading modifyMeshDict\n" << endl;
00366 
00367     IOdictionary dict
00368     (
00369         IOobject
00370         (
00371             "modifyMeshDict",
00372             runTime.system(),
00373             mesh,
00374             IOobject::MUST_READ,
00375             IOobject::NO_WRITE
00376         )
00377     );
00378 
00379     // Read all from the dictionary.
00380     List<Pair<point> > pointsToMove(dict.lookup("pointsToMove"));
00381     List<Pair<point> > edgesToSplit(dict.lookup("edgesToSplit"));
00382     List<Pair<point> > facesToTriangulate
00383     (
00384         dict.lookup("facesToTriangulate")
00385     );
00386 
00387     bool cutBoundary =
00388     (
00389         pointsToMove.size()
00390      || edgesToSplit.size()
00391      || facesToTriangulate.size()
00392     );
00393 
00394     List<Pair<point> > edgesToCollapse(dict.lookup("edgesToCollapse"));
00395 
00396     bool collapseEdge = edgesToCollapse.size();
00397 
00398     List<Pair<point> > cellsToPyramidise(dict.lookup("cellsToSplit"));
00399 
00400     bool cellsToSplit = cellsToPyramidise.size();
00401 
00402     //List<Tuple<pointField,point> >
00403     //  cellsToCreate(dict.lookup("cellsToCreate"));
00404 
00405     Info<< "Read from " << dict.name() << nl
00406         << "  Boundary cutting module:" << nl
00407         << "    points to move      :" << pointsToMove.size() << nl
00408         << "    edges to split      :" << edgesToSplit.size() << nl
00409         << "    faces to triangulate:" << facesToTriangulate.size() << nl
00410         << "  Cell splitting module:" << nl
00411         << "    cells to split      :" << cellsToPyramidise.size() << nl
00412         << "  Edge collapsing module:" << nl
00413         << "    edges to collapse   :" << edgesToCollapse.size() << nl
00414         //<< "    cells to create     :" << cellsToCreate.size() << nl
00415         << endl;
00416 
00417     if
00418     (
00419         (cutBoundary && collapseEdge)
00420      || (cutBoundary && cellsToSplit)
00421      || (collapseEdge && cellsToSplit)
00422     )
00423     {
00424         FatalErrorIn(args.executable())
00425             << "Used more than one mesh modifying module "
00426             << "(boundary cutting, cell splitting, edge collapsing)" << nl
00427             << "Please do them in separate passes." << exit(FatalError);
00428     }
00429 
00430 
00431 
00432     // Get calculating engine for all of outside
00433     const SubList<face> outsideFaces
00434     (
00435         mesh.faces(),
00436         mesh.nFaces() - mesh.nInternalFaces(),
00437         mesh.nInternalFaces()
00438     );
00439 
00440     primitivePatch allBoundary(outsideFaces, mesh.points());
00441 
00442 
00443     // Look up mesh labels and convert to input for boundaryCutter.
00444 
00445     bool validInputs = true;
00446 
00447 
00448     Info<< nl << "Looking up points to move ..." << nl << endl;
00449     Map<point> pointToPos(pointsToMove.size());
00450     forAll(pointsToMove, i)
00451     {
00452         const Pair<point>& pts = pointsToMove[i];
00453 
00454         label pointI = findPoint(allBoundary, pts.first());
00455 
00456         if (pointI == -1 || !pointToPos.insert(pointI, pts.second()))
00457         {
00458             Info<< "Could not insert mesh point " << pointI
00459                 << " for input point " << pts.first() << nl
00460                 << "Perhaps the point is already marked for moving?" << endl;
00461             validInputs = false;
00462         }
00463     }
00464 
00465 
00466     Info<< nl << "Looking up edges to split ..." << nl << endl;
00467     Map<List<point> > edgeToCuts(edgesToSplit.size());
00468     forAll(edgesToSplit, i)
00469     {
00470         const Pair<point>& pts = edgesToSplit[i];
00471 
00472         label edgeI = findEdge(mesh, allBoundary, pts.first());
00473 
00474         if
00475         (
00476             edgeI == -1
00477         || !edgeToCuts.insert(edgeI, List<point>(1, pts.second()))
00478         )
00479         {
00480             Info<< "Could not insert mesh edge " << edgeI
00481                 << " for input point " << pts.first() << nl
00482                 << "Perhaps the edge is already marked for cutting?" << endl;
00483 
00484             validInputs = false;
00485         }
00486     }
00487 
00488 
00489     Info<< nl << "Looking up faces to triangulate ..." << nl << endl;
00490     Map<point> faceToDecompose(facesToTriangulate.size());
00491     forAll(facesToTriangulate, i)
00492     {
00493         const Pair<point>& pts = facesToTriangulate[i];
00494 
00495         label faceI = findFace(mesh, allBoundary, pts.first());
00496 
00497         if (faceI == -1 || !faceToDecompose.insert(faceI, pts.second()))
00498         {
00499             Info<< "Could not insert mesh face " << faceI
00500                 << " for input point " << pts.first() << nl
00501                 << "Perhaps the face is already marked for splitting?" << endl;
00502 
00503             validInputs = false;
00504         }
00505     }
00506 
00507 
00508 
00509     Info<< nl << "Looking up cells to convert to pyramids around"
00510         << " cell centre ..." << nl << endl;
00511     Map<point> cellToPyrCentre(cellsToPyramidise.size());
00512     forAll(cellsToPyramidise, i)
00513     {
00514         const Pair<point>& pts = cellsToPyramidise[i];
00515 
00516         label cellI = findCell(mesh, pts.first());
00517 
00518         if (cellI == -1 || !cellToPyrCentre.insert(cellI, pts.second()))
00519         {
00520             Info<< "Could not insert mesh cell " << cellI
00521                 << " for input point " << pts.first() << nl
00522                 << "Perhaps the cell is already marked for splitting?" << endl;
00523 
00524             validInputs = false;
00525         }
00526     }
00527 
00528 
00529     Info<< nl << "Looking up edges to collapse ..." << nl << endl;
00530     Map<point> edgeToPos(edgesToCollapse.size());
00531     forAll(edgesToCollapse, i)
00532     {
00533         const Pair<point>& pts = edgesToCollapse[i];
00534 
00535         label edgeI = findEdge(mesh, allBoundary, pts.first());
00536 
00537         if (edgeI == -1 || !edgeToPos.insert(edgeI, pts.second()))
00538         {
00539             Info<< "Could not insert mesh edge " << edgeI
00540                 << " for input point " << pts.first() << nl
00541                 << "Perhaps the edge is already marked for collaping?" << endl;
00542 
00543             validInputs = false;
00544         }
00545     }
00546 
00547 
00548 
00549     if (!validInputs)
00550     {
00551         Info<< nl << "There was a problem in one of the inputs in the"
00552             << " dictionary. Not modifying mesh." << endl;
00553     }
00554     else if (cellToPyrCentre.size())
00555     {
00556         Info<< nl << "All input cells located. Modifying mesh." << endl;
00557 
00558         // Mesh change engine
00559         cellSplitter cutter(mesh);
00560 
00561         // Topo change container
00562         polyTopoChange meshMod(mesh);
00563 
00564         // Insert commands into meshMod
00565         cutter.setRefinement(cellToPyrCentre, meshMod);
00566 
00567         // Do changes
00568         autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
00569 
00570         if (morphMap().hasMotionPoints())
00571         {
00572             mesh.movePoints(morphMap().preMotionPoints());
00573         }
00574 
00575         cutter.updateMesh(morphMap());
00576 
00577         if (!overwrite)
00578         {
00579             runTime++;
00580         }
00581         else
00582         {
00583             mesh.setInstance(oldInstance);
00584         }
00585 
00586         // Write resulting mesh
00587         Info << "Writing modified mesh to time " << runTime.timeName() << endl;
00588         mesh.write();
00589     }
00590     else if (edgeToPos.size())
00591     {
00592         Info<< nl << "All input edges located. Modifying mesh." << endl;
00593 
00594         // Mesh change engine
00595         edgeCollapser cutter(mesh);
00596 
00597         pointField newPoints(mesh.points());
00598 
00599         // Get new positions and construct collapse network
00600         forAllConstIter(Map<point>, edgeToPos, iter)
00601         {
00602             label edgeI = iter.key();
00603             const edge& e = mesh.edges()[edgeI];
00604 
00605             cutter.collapseEdge(edgeI, e[0]);
00606             newPoints[e[0]] = iter();
00607         }
00608 
00609         // Move master point to destination.
00610         mesh.movePoints(newPoints);
00611 
00612         // Topo change container
00613         polyTopoChange meshMod(mesh);
00614 
00615         // Insert
00616         cutter.setRefinement(meshMod);
00617 
00618         // Do changes
00619         autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
00620 
00621         if (morphMap().hasMotionPoints())
00622         {
00623             mesh.movePoints(morphMap().preMotionPoints());
00624         }
00625 
00626         // Not implemented yet:
00627         //cutter.updateMesh(morphMap());
00628 
00629 
00630         if (!overwrite)
00631         {
00632             runTime++;
00633         }
00634         else
00635         {
00636             mesh.setInstance(oldInstance);
00637         }
00638 
00639         // Write resulting mesh
00640         Info << "Writing modified mesh to time " << runTime.timeName() << endl;
00641         mesh.write();
00642     }
00643     else
00644     {
00645         Info<< nl << "All input points located. Modifying mesh." << endl;
00646 
00647         // Mesh change engine
00648         boundaryCutter cutter(mesh);
00649 
00650         // Topo change container
00651         polyTopoChange meshMod(mesh);
00652 
00653         // Insert commands into meshMod
00654         cutter.setRefinement
00655         (
00656             pointToPos,
00657             edgeToCuts,
00658             Map<labelPair>(0),  // Faces to split diagonally
00659             faceToDecompose,    // Faces to triangulate
00660             meshMod
00661         );
00662 
00663         // Do changes
00664         autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
00665 
00666         if (morphMap().hasMotionPoints())
00667         {
00668             mesh.movePoints(morphMap().preMotionPoints());
00669         }
00670 
00671         cutter.updateMesh(morphMap());
00672 
00673         if (!overwrite)
00674         {
00675             runTime++;
00676         }
00677         else
00678         {
00679             mesh.setInstance(oldInstance);
00680         }
00681 
00682         // Write resulting mesh
00683         Info << "Writing modified mesh to time " << runTime.timeName() << endl;
00684         mesh.write();
00685     }
00686 
00687 
00688     Info << nl << "End" << endl;
00689 
00690     return 0;
00691 }
00692 
00693 
00694 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines