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

extrudedMesh.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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "extrudedMesh.H"
00027 #include <OpenFOAM/wallPolyPatch.H>
00028 #include <meshTools/meshTools.H>
00029 #include <OpenFOAM/ListOps.H>
00030 #include <OpenFOAM/OFstream.H>
00031 
00032 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00033 
00034 bool Foam::extrudedMesh::sameOrder(const face& f, const edge& e)
00035 {
00036     label i = findIndex(f, e[0]);
00037 
00038     label nextI = (i == f.size()-1 ? 0 : i+1);
00039 
00040     return f[nextI] == e[1];
00041 }
00042 
00043 
00044 template
00045 <
00046     class Face,
00047     template<class> class FaceList,
00048     class PointField
00049 >
00050 Foam::Xfer<Foam::pointField> Foam::extrudedMesh::extrudedPoints
00051 (
00052     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
00053     const extrudeModel& model
00054 )
00055 {
00056     const pointField& surfacePoints = extrudePatch.localPoints();
00057     const vectorField& surfaceNormals = extrudePatch.pointNormals();
00058 
00059     const label nLayers = model.nLayers();
00060 
00061     pointField ePoints((nLayers + 1)*surfacePoints.size());
00062 
00063     for (label layer=0; layer<=nLayers; layer++)
00064     {
00065         label offset = layer*surfacePoints.size();
00066 
00067         forAll(surfacePoints, i)
00068         {
00069             ePoints[offset + i] = model
00070             (
00071                 surfacePoints[i],
00072                 surfaceNormals[i],
00073                 layer
00074             );
00075         }
00076     }
00077 
00078     // return points for transferring
00079     return xferMove(ePoints);
00080 }
00081 
00082 
00083 template<class Face, template<class> class FaceList, class PointField>
00084 Foam::Xfer<Foam::faceList> Foam::extrudedMesh::extrudedFaces
00085 (
00086     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
00087     const extrudeModel& model
00088 )
00089 {
00090     const pointField& surfacePoints = extrudePatch.localPoints();
00091     const List<face>& surfaceFaces = extrudePatch.localFaces();
00092     const edgeList& surfaceEdges = extrudePatch.edges();
00093     const label nInternalEdges = extrudePatch.nInternalEdges();
00094 
00095     const label nLayers = model.nLayers();
00096 
00097     label nFaces = 
00098         (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
00099 
00100     faceList eFaces(nFaces);
00101 
00102     labelList quad(4);
00103     label facei = 0;
00104 
00105     // Internal faces
00106     for (label layer=0; layer<nLayers; layer++)
00107     {
00108         label currentLayerOffset = layer*surfacePoints.size();
00109         label nextLayerOffset = currentLayerOffset + surfacePoints.size();
00110 
00111         // Vertical faces from layer to layer+1
00112         for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
00113         {
00114             const edge& e = surfaceEdges[edgeI];
00115             const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
00116 
00117             face& f = eFaces[facei++];
00118             f.setSize(4);
00119 
00120             if
00121             (
00122                 (edgeFaces[0] < edgeFaces[1])
00123              == sameOrder(surfaceFaces[edgeFaces[0]], e)
00124             )
00125             {
00126                 f[0] = e[0] + currentLayerOffset;
00127                 f[1] = e[1] + currentLayerOffset;
00128                 f[2] = e[1] + nextLayerOffset;
00129                 f[3] = e[0] + nextLayerOffset;
00130             }
00131             else
00132             {
00133                 f[0] = e[1] + currentLayerOffset;
00134                 f[1] = e[0] + currentLayerOffset;
00135                 f[2] = e[0] + nextLayerOffset;
00136                 f[3] = e[1] + nextLayerOffset;
00137             }
00138         }
00139 
00140         // Faces between layer and layer+1
00141         if (layer < nLayers-1)
00142         {
00143             forAll(surfaceFaces, i)
00144             {
00145                 eFaces[facei++] =
00146                     face
00147                     (
00148                         surfaceFaces[i] //.reverseFace()
00149                       + nextLayerOffset
00150                     );
00151             }
00152         }
00153     }
00154 
00155     // External side faces
00156     for (label layer=0; layer<nLayers; layer++)
00157     {
00158         label currentLayerOffset = layer*surfacePoints.size();
00159         label nextLayerOffset = currentLayerOffset + surfacePoints.size();
00160 
00161         // Side faces across layer
00162         for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
00163         {
00164             const edge& e = surfaceEdges[edgeI];
00165             const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
00166 
00167             face& f = eFaces[facei++];
00168             f.setSize(4);
00169 
00170             if (sameOrder(surfaceFaces[edgeFaces[0]], e))
00171             {
00172                 f[0] = e[0] + currentLayerOffset;
00173                 f[1] = e[1] + currentLayerOffset;
00174                 f[2] = e[1] + nextLayerOffset;
00175                 f[3] = e[0] + nextLayerOffset;
00176             }
00177             else
00178             {
00179                 f[0] = e[1] + currentLayerOffset;
00180                 f[1] = e[0] + currentLayerOffset;
00181                 f[2] = e[0] + nextLayerOffset;
00182                 f[3] = e[1] + nextLayerOffset;
00183             }
00184         }
00185     }
00186 
00187     // Bottom faces
00188     forAll(surfaceFaces, i)
00189     {
00190         eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
00191     }
00192 
00193     // Top faces
00194     forAll(surfaceFaces, i)
00195     {
00196         eFaces[facei++] =
00197             face
00198             (
00199                 surfaceFaces[i]
00200               + nLayers*surfacePoints.size()
00201             );
00202     }
00203 
00204     // return points for transferring
00205     return xferMove(eFaces);
00206 }
00207 
00208 
00209 template<class Face, template<class> class FaceList, class PointField>
00210 Foam::Xfer<Foam::cellList> Foam::extrudedMesh::extrudedCells
00211 (
00212     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
00213     const extrudeModel& model
00214 )
00215 {
00216     const List<face>& surfaceFaces = extrudePatch.localFaces();
00217     const edgeList& surfaceEdges = extrudePatch.edges();
00218     const label nInternalEdges = extrudePatch.nInternalEdges();
00219 
00220     const label nLayers = model.nLayers();
00221 
00222     cellList eCells(nLayers*surfaceFaces.size());
00223 
00224     // Size the cells
00225     forAll(surfaceFaces, i)
00226     {
00227         const face& f = surfaceFaces[i];
00228 
00229         for (label layer=0; layer<nLayers; layer++)
00230         {
00231             eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
00232         }
00233     }
00234 
00235     // Current face count per cell.
00236     labelList nCellFaces(eCells.size(), 0);
00237 
00238 
00239     label facei = 0;
00240 
00241     for (label layer=0; layer<nLayers; layer++)
00242     {
00243         // Side faces from layer to layer+1
00244         for (label i=0; i<nInternalEdges; i++)
00245         {
00246             // Get patch faces using edge
00247             const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
00248 
00249             // Get cells on this layer
00250             label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
00251             label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
00252 
00253             eCells[cell0][nCellFaces[cell0]++] = facei;
00254             eCells[cell1][nCellFaces[cell1]++] = facei;
00255 
00256             facei++;
00257         }
00258 
00259         // Faces between layer and layer+1
00260         if (layer < nLayers-1)
00261         {
00262             forAll(surfaceFaces, i)
00263             {
00264                 label cell0 = layer*surfaceFaces.size() + i;
00265                 label cell1 = (layer+1)*surfaceFaces.size() + i;
00266 
00267                 eCells[cell0][nCellFaces[cell0]++] = facei;
00268                 eCells[cell1][nCellFaces[cell1]++] = facei;
00269 
00270                 facei++;
00271             }
00272         }
00273     }
00274 
00275     // External side faces
00276     for (label layer=0; layer<nLayers; layer++)
00277     {
00278         // Side faces across layer
00279         for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
00280         {
00281             // Get patch faces using edge
00282             const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
00283 
00284             // Get cells on this layer
00285             label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
00286 
00287             eCells[cell0][nCellFaces[cell0]++] = facei;
00288 
00289             facei++;
00290         }
00291     }
00292 
00293     // Top faces
00294     forAll(surfaceFaces, i)
00295     {
00296         eCells[i][nCellFaces[i]++] = facei;
00297 
00298         facei++;
00299     }
00300 
00301     // Bottom faces
00302     forAll(surfaceFaces, i)
00303     {
00304         label cell0 = (nLayers-1)*surfaceFaces.size() + i;
00305 
00306         eCells[cell0][nCellFaces[cell0]++] = facei;
00307 
00308         facei++;
00309     }
00310 
00311     // return points for transferring
00312     return xferMove(eCells);
00313 }
00314 
00315 
00316 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00317 
00318 template
00319 <
00320     class Face,
00321     template<class> class FaceList,
00322     class PointField
00323 >
00324 Foam::extrudedMesh::extrudedMesh
00325 (
00326     const IOobject& io,
00327     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
00328     const extrudeModel& model
00329 )
00330 :
00331     polyMesh
00332     (
00333         io,
00334         extrudedPoints(extrudePatch, model),
00335         extrudedFaces(extrudePatch, model),
00336         extrudedCells(extrudePatch, model)
00337     ),
00338     model_(model)
00339 {
00340     List<polyPatch*> patches(3);
00341 
00342     label facei = nInternalFaces();
00343 
00344     label sz =
00345         model_.nLayers()
00346        *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
00347 
00348     patches[0] = new wallPolyPatch
00349     (
00350         "sides",
00351         sz,
00352         facei,
00353         0,
00354         boundaryMesh()
00355     );
00356 
00357     facei += sz;
00358 
00359     patches[1] = new polyPatch
00360     (
00361         "originalPatch",
00362         extrudePatch.size(),
00363         facei,
00364         1,
00365         boundaryMesh()
00366     );
00367 
00368     facei += extrudePatch.size();
00369 
00370     patches[2] = new polyPatch
00371     (
00372         "otherSide",
00373         extrudePatch.size(),
00374         facei,
00375         2,
00376         boundaryMesh()
00377     );
00378 
00379     addPatches(patches);
00380 }
00381 
00382 
00383 // ************************ vim: set sw=4 sts=4 et: ************************ //
00384 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines