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

refineHexMesh.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     refineHexMesh
00026 
00027 Description
00028     Refines a hex mesh by 2x2x2 cell splitting.
00029 
00030 Usage
00031 
00032     - refineHexMesh [OPTIONS] <cellSet name>
00033 
00034     @param <cellSet name> \n
00035     @todo Detailed description of argument.
00036 
00037     @param -overwrite \n
00038     Overwrite existing data.
00039 
00040     @param -case <dir>\n
00041     Case directory.
00042 
00043     @param -parallel \n
00044     Run in parallel.
00045 
00046     @param -help \n
00047     Display help message.
00048 
00049     @param -doc \n
00050     Display Doxygen API documentation page for this application.
00051 
00052     @param -srcDoc \n
00053     Display Doxygen source documentation page for this application.
00054 
00055 \*---------------------------------------------------------------------------*/
00056 
00057 #include <finiteVolume/fvMesh.H>
00058 #include <OpenFOAM/pointMesh.H>
00059 #include <OpenFOAM/argList.H>
00060 #include <OpenFOAM/Time.H>
00061 #include <dynamicMesh/hexRef8.H>
00062 #include <meshTools/cellSet.H>
00063 #include <OpenFOAM/OFstream.H>
00064 #include <meshTools/meshTools.H>
00065 #include <OpenFOAM/IFstream.H>
00066 #include <dynamicMesh/polyTopoChange.H>
00067 #include <OpenFOAM/mapPolyMesh.H>
00068 #include <finiteVolume/volMesh.H>
00069 #include <finiteVolume/surfaceMesh.H>
00070 #include <finiteVolume/volFields.H>
00071 #include <finiteVolume/surfaceFields.H>
00072 #include <OpenFOAM/pointFields.H>
00073 #include <OpenFOAM/ReadFields.H>
00074 
00075 using namespace Foam;
00076 
00077 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00078 
00079 // Main program:
00080 int main(int argc, char *argv[])
00081 {
00082     argList::validOptions.insert("overwrite", "");
00083     argList::validArgs.append("cellSet");
00084 #   include <OpenFOAM/setRootCase.H>
00085 #   include <OpenFOAM/createTime.H>
00086     runTime.functionObjects().off();
00087 #   include <OpenFOAM/createMesh.H>
00088     const word oldInstance = mesh.pointsInstance();
00089 
00090     pointMesh pMesh(mesh);
00091 
00092     word cellSetName(args.args()[1]);
00093     bool overwrite = args.optionFound("overwrite");
00094 
00095     Info<< "Reading cells to refine from cellSet " << cellSetName
00096         << nl << endl;
00097 
00098     cellSet cellsToRefine(mesh, cellSetName);
00099 
00100     Info<< "Read " << returnReduce(cellsToRefine.size(), sumOp<label>())
00101         << " cells to refine from cellSet " << cellSetName << nl
00102         << endl;
00103 
00104 
00105     // Read objects in time directory
00106     IOobjectList objects(mesh, runTime.timeName());
00107 
00108     // Read vol fields.
00109 
00110     PtrList<volScalarField> vsFlds;
00111     ReadFields(mesh, objects, vsFlds);
00112 
00113     PtrList<volVectorField> vvFlds;
00114     ReadFields(mesh, objects, vvFlds);
00115 
00116     PtrList<volSphericalTensorField> vstFlds;
00117     ReadFields(mesh, objects, vstFlds);
00118 
00119     PtrList<volSymmTensorField> vsymtFlds;
00120     ReadFields(mesh, objects, vsymtFlds);
00121 
00122     PtrList<volTensorField> vtFlds;
00123     ReadFields(mesh, objects, vtFlds);
00124 
00125     // Read surface fields.
00126 
00127     PtrList<surfaceScalarField> ssFlds;
00128     ReadFields(mesh, objects, ssFlds);
00129 
00130     PtrList<surfaceVectorField> svFlds;
00131     ReadFields(mesh, objects, svFlds);
00132 
00133     PtrList<surfaceSphericalTensorField> sstFlds;
00134     ReadFields(mesh, objects, sstFlds);
00135 
00136     PtrList<surfaceSymmTensorField> ssymtFlds;
00137     ReadFields(mesh, objects, ssymtFlds);
00138 
00139     PtrList<surfaceTensorField> stFlds;
00140     ReadFields(mesh, objects, stFlds);
00141 
00142     // Read point fields
00143     PtrList<pointScalarField> psFlds;
00144     ReadFields(pMesh, objects, psFlds);
00145 
00146     PtrList<pointVectorField> pvFlds;
00147     ReadFields(pMesh, objects, pvFlds);
00148 
00149 
00150 
00151     // Construct refiner without unrefinement. Read existing point/cell level.
00152     hexRef8 meshCutter(mesh);
00153 
00154     // Some stats
00155     Info<< "Read mesh:" << nl
00156         << "    cells:" << mesh.globalData().nTotalCells() << nl
00157         << "    faces:" << mesh.globalData().nTotalFaces() << nl
00158         << "    points:" << mesh.globalData().nTotalPoints() << nl
00159         << "    cellLevel :"
00160         << " min:" << gMin(meshCutter.cellLevel())
00161         << " max:" << gMax(meshCutter.cellLevel()) << nl
00162         << "    pointLevel :"
00163         << " min:" << gMin(meshCutter.pointLevel())
00164         << " max:" << gMax(meshCutter.pointLevel()) << nl
00165         << endl;
00166 
00167 
00168     // Maintain 2:1 ratio
00169     labelList newCellsToRefine
00170     (
00171         meshCutter.consistentRefinement
00172         (
00173             cellsToRefine.toc(),
00174             true                  // extend set
00175         )
00176     );
00177 
00178     // Mesh changing engine.
00179     polyTopoChange meshMod(mesh);
00180 
00181     // Play refinement commands into mesh changer.
00182     meshCutter.setRefinement(newCellsToRefine, meshMod);
00183 
00184     if (!overwrite)
00185     {
00186         runTime++;
00187     }
00188 
00189     // Create mesh, return map from old to new mesh.
00190     autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
00191 
00192     // Update fields
00193     mesh.updateMesh(map);
00194     pMesh.updateMesh(map);
00195 
00196     // Update numbering of cells/vertices.
00197     meshCutter.updateMesh(map);
00198 
00199     // Optionally inflate mesh
00200     if (map().hasMotionPoints())
00201     {
00202         mesh.movePoints(map().preMotionPoints());
00203         pMesh.movePoints(map().preMotionPoints());
00204     }
00205 
00206     Pout<< "Refined from " << returnReduce(map().nOldCells(), sumOp<label>())
00207         << " to " << mesh.globalData().nTotalCells() << " cells." << nl << endl;
00208 
00209     if (overwrite)
00210     {
00211         mesh.setInstance(oldInstance);
00212         meshCutter.setInstance(oldInstance);
00213     }
00214     Info<< "Writing mesh to " << runTime.timeName() << endl;
00215 
00216     mesh.write();
00217     meshCutter.write();
00218 
00219     Info<< "End\n" << endl;
00220 
00221     return 0;
00222 }
00223 
00224 
00225 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines