fix face direction
[OpenFOAM-1.5.x.git] / applications / utilities / mesh / generation / extrudeMesh / extrudedMesh / extrudedMesh.C
blob1e84cfa689e4e6daf439fe14d3b05e4bb2f32920
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "extrudedMesh.H"
28 #include "wallPolyPatch.H"
29 #include "meshTools.H"
30 #include "ListOps.H"
31 #include "OFstream.H"
33 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
35 bool Foam::extrudedMesh::sameOrder(const face& f, const edge& e)
37     label i = findIndex(f, e[0]);
39     label nextI = (i == f.size()-1 ? 0 : i+1);
41     return f[nextI] == e[1];
45 template
47     class Face,
48     template<class> class FaceList,
49     class PointField
51 Foam::pointField Foam::extrudedMesh::extrudedPoints
53     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
54     const extrudeModel& model
57     const pointField& surfacePoints = extrudePatch.localPoints();
58     const vectorField& surfaceNormals = extrudePatch.pointNormals();
60     const label nLayers = model.nLayers();
62     pointField ePoints((nLayers + 1)*surfacePoints.size());
64     for (label layer=0; layer<=nLayers; layer++)
65     {
66         label offset = layer*surfacePoints.size();
68         forAll(surfacePoints, i)
69         {
70             ePoints[offset + i] = model
71             (
72                 surfacePoints[i],
73                 surfaceNormals[i],
74                 layer
75             );
76         }
77     }
79     return ePoints;
83 template<class Face, template<class> class FaceList, class PointField>
84 Foam::faceList Foam::extrudedMesh::extrudedFaces
86     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
87     const extrudeModel& model
90     const pointField& surfacePoints = extrudePatch.localPoints();
91     const List<face>& surfaceFaces = extrudePatch.localFaces();
92     const edgeList& surfaceEdges = extrudePatch.edges();
93     const label nInternalEdges = extrudePatch.nInternalEdges();
95     const label nLayers = model.nLayers();
97     label nFaces = 
98         (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
100     faceList eFaces(nFaces);
102     labelList quad(4);
103     label facei = 0;
105     // Internal faces
106     for (label layer=0; layer<nLayers; layer++)
107     {
108         label currentLayerOffset = layer*surfacePoints.size();
109         label nextLayerOffset = currentLayerOffset + surfacePoints.size();
111         // Vertical faces from layer to layer+1
112         for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
113         {
114             const edge& e = surfaceEdges[edgeI];
115             const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
117             face& f = eFaces[facei++];
118             f.setSize(4);
120             if
121             (
122                 (edgeFaces[0] < edgeFaces[1])
123              == sameOrder(surfaceFaces[edgeFaces[0]], e)
124             )
125             {
126                 f[0] = e[0] + currentLayerOffset;
127                 f[1] = e[1] + currentLayerOffset;
128                 f[2] = e[1] + nextLayerOffset;
129                 f[3] = e[0] + nextLayerOffset;
130             }
131             else
132             {
133                 f[0] = e[1] + currentLayerOffset;
134                 f[1] = e[0] + currentLayerOffset;
135                 f[2] = e[0] + nextLayerOffset;
136                 f[3] = e[1] + nextLayerOffset;
137             }
138         }
140         // Faces between layer and layer+1
141         if (layer < nLayers-1)
142         {
143             forAll(surfaceFaces, i)
144             {
145                 eFaces[facei++] =
146                     face
147                     (
148                         surfaceFaces[i] //.reverseFace()
149                       + nextLayerOffset
150                     );
151             }
152         }
153     }
155     // External side faces
156     for (label layer=0; layer<nLayers; layer++)
157     {
158         label currentLayerOffset = layer*surfacePoints.size();
159         label nextLayerOffset = currentLayerOffset + surfacePoints.size();
161         // Side faces across layer
162         for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
163         {
164             const edge& e = surfaceEdges[edgeI];
165             const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
167             face& f = eFaces[facei++];
168             f.setSize(4);
170             if (sameOrder(surfaceFaces[edgeFaces[0]], e))
171             {
172                 f[0] = e[0] + currentLayerOffset;
173                 f[1] = e[1] + currentLayerOffset;
174                 f[2] = e[1] + nextLayerOffset;
175                 f[3] = e[0] + nextLayerOffset;
176             }
177             else
178             {
179                 f[0] = e[1] + currentLayerOffset;
180                 f[1] = e[0] + currentLayerOffset;
181                 f[2] = e[0] + nextLayerOffset;
182                 f[3] = e[1] + nextLayerOffset;
183             }
184         }
185     }
187     // Bottom faces
188     forAll(surfaceFaces, i)
189     {
190         eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
191     }
193     // Top faces
194     forAll(surfaceFaces, i)
195     {
196         eFaces[facei++] =
197             face
198             (
199                 surfaceFaces[i]
200               + nLayers*surfacePoints.size()
201             );
202     }
204     return eFaces;
208 template<class Face, template<class> class FaceList, class PointField>
209 Foam::cellList Foam::extrudedMesh::extrudedCells
211     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
212     const extrudeModel& model
215     const List<face>& surfaceFaces = extrudePatch.localFaces();
216     const edgeList& surfaceEdges = extrudePatch.edges();
217     const label nInternalEdges = extrudePatch.nInternalEdges();
219     const label nLayers = model.nLayers();
221     cellList eCells(nLayers*surfaceFaces.size());
223     // Size the cells
224     forAll(surfaceFaces, i)
225     {
226         const face& f = surfaceFaces[i];
228         for (label layer=0; layer<nLayers; layer++)
229         {
230             eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
231         }
232     }
234     // Current face count per cell.
235     labelList nCellFaces(eCells.size(), 0);
238     label facei = 0;
240     for (label layer=0; layer<nLayers; layer++)
241     {
242         // Side faces from layer to layer+1
243         for (label i=0; i<nInternalEdges; i++)
244         {
245             // Get patch faces using edge
246             const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
248             // Get cells on this layer
249             label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
250             label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
252             eCells[cell0][nCellFaces[cell0]++] = facei;
253             eCells[cell1][nCellFaces[cell1]++] = facei;
255             facei++;
256         }
258         // Faces between layer and layer+1
259         if (layer < nLayers-1)
260         {
261             forAll(surfaceFaces, i)
262             {
263                 label cell0 = layer*surfaceFaces.size() + i;
264                 label cell1 = (layer+1)*surfaceFaces.size() + i;
266                 eCells[cell0][nCellFaces[cell0]++] = facei;
267                 eCells[cell1][nCellFaces[cell1]++] = facei;
269                 facei++;
270             }
271         }
272     }
274     // External side faces
275     for (label layer=0; layer<nLayers; layer++)
276     {
277         // Side faces across layer
278         for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
279         {
280             // Get patch faces using edge
281             const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
283             // Get cells on this layer
284             label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
286             eCells[cell0][nCellFaces[cell0]++] = facei;
288             facei++;
289         }
290     }
292     // Top faces
293     forAll(surfaceFaces, i)
294     {
295         eCells[i][nCellFaces[i]++] = facei;
297         facei++;
298     }
300     // Bottom faces
301     forAll(surfaceFaces, i)
302     {
303         label cell0 = (nLayers-1)*surfaceFaces.size() + i;
305         eCells[cell0][nCellFaces[cell0]++] = facei;
307         facei++;
308     }
310     return eCells;
314 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
316 template
318     class Face,
319     template<class> class FaceList,
320     class PointField
322 Foam::extrudedMesh::extrudedMesh
324     const IOobject& io,
325     const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
326     const extrudeModel& model
329     polyMesh
330     (
331         io,
332         extrudedPoints(extrudePatch, model),
333         extrudedFaces(extrudePatch, model),
334         extrudedCells(extrudePatch, model)
335     ),
336     model_(model)
338     List<polyPatch*> patches(3);
340     label facei = nInternalFaces();
342     label sz =
343         model_.nLayers()
344        *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
346     patches[0] = new wallPolyPatch
347     (
348         "sides",
349         sz,
350         facei,
351         0,
352         boundaryMesh()
353     );
355     facei += sz;
357     patches[1] = new polyPatch
358     (
359         "originalPatch",
360         extrudePatch.size(),
361         facei,
362         1,
363         boundaryMesh()
364     );
366     facei += extrudePatch.size();
368     patches[2] = new polyPatch
369     (
370         "otherSide",
371         extrudePatch.size(),
372         facei,
373         2,
374         boundaryMesh()
375     );
377     addPatches(patches);
381 // ************************************************************************* //