1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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"
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];
48 template<class> class FaceList,
51 Foam::Xfer<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++)
66 label offset = layer*surfacePoints.size();
68 forAll(surfacePoints, i)
70 ePoints[offset + i] = model
79 // return points for transferring
80 return xferMove(ePoints);
84 template<class Face, template<class> class FaceList, class PointField>
85 Foam::Xfer<Foam::faceList> Foam::extrudedMesh::extrudedFaces
87 const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
88 const extrudeModel& model
91 const pointField& surfacePoints = extrudePatch.localPoints();
92 const List<face>& surfaceFaces = extrudePatch.localFaces();
93 const edgeList& surfaceEdges = extrudePatch.edges();
94 const label nInternalEdges = extrudePatch.nInternalEdges();
96 const label nLayers = model.nLayers();
99 (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
101 faceList eFaces(nFaces);
107 for (label layer=0; layer<nLayers; layer++)
109 label currentLayerOffset = layer*surfacePoints.size();
110 label nextLayerOffset = currentLayerOffset + surfacePoints.size();
112 // Vertical faces from layer to layer+1
113 for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
115 const edge& e = surfaceEdges[edgeI];
116 const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
118 face& f = eFaces[facei++];
123 (edgeFaces[0] < edgeFaces[1])
124 == sameOrder(surfaceFaces[edgeFaces[0]], e)
127 f[0] = e[0] + currentLayerOffset;
128 f[1] = e[1] + currentLayerOffset;
129 f[2] = e[1] + nextLayerOffset;
130 f[3] = e[0] + nextLayerOffset;
134 f[0] = e[1] + currentLayerOffset;
135 f[1] = e[0] + currentLayerOffset;
136 f[2] = e[0] + nextLayerOffset;
137 f[3] = e[1] + nextLayerOffset;
141 // Faces between layer and layer+1
142 if (layer < nLayers-1)
144 forAll(surfaceFaces, i)
149 surfaceFaces[i] //.reverseFace()
156 // External side faces
157 for (label layer=0; layer<nLayers; layer++)
159 label currentLayerOffset = layer*surfacePoints.size();
160 label nextLayerOffset = currentLayerOffset + surfacePoints.size();
162 // Side faces across layer
163 for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
165 const edge& e = surfaceEdges[edgeI];
166 const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
168 face& f = eFaces[facei++];
171 if (sameOrder(surfaceFaces[edgeFaces[0]], e))
173 f[0] = e[0] + currentLayerOffset;
174 f[1] = e[1] + currentLayerOffset;
175 f[2] = e[1] + nextLayerOffset;
176 f[3] = e[0] + nextLayerOffset;
180 f[0] = e[1] + currentLayerOffset;
181 f[1] = e[0] + currentLayerOffset;
182 f[2] = e[0] + nextLayerOffset;
183 f[3] = e[1] + nextLayerOffset;
189 forAll(surfaceFaces, i)
191 eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
195 forAll(surfaceFaces, i)
201 + nLayers*surfacePoints.size()
205 // return points for transferring
206 return xferMove(eFaces);
210 template<class Face, template<class> class FaceList, class PointField>
211 Foam::Xfer<Foam::cellList> Foam::extrudedMesh::extrudedCells
213 const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
214 const extrudeModel& model
217 const List<face>& surfaceFaces = extrudePatch.localFaces();
218 const edgeList& surfaceEdges = extrudePatch.edges();
219 const label nInternalEdges = extrudePatch.nInternalEdges();
221 const label nLayers = model.nLayers();
223 cellList eCells(nLayers*surfaceFaces.size());
226 forAll(surfaceFaces, i)
228 const face& f = surfaceFaces[i];
230 for (label layer=0; layer<nLayers; layer++)
232 eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
236 // Current face count per cell.
237 labelList nCellFaces(eCells.size(), 0);
242 for (label layer=0; layer<nLayers; layer++)
244 // Side faces from layer to layer+1
245 for (label i=0; i<nInternalEdges; i++)
247 // Get patch faces using edge
248 const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
250 // Get cells on this layer
251 label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
252 label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
254 eCells[cell0][nCellFaces[cell0]++] = facei;
255 eCells[cell1][nCellFaces[cell1]++] = facei;
260 // Faces between layer and layer+1
261 if (layer < nLayers-1)
263 forAll(surfaceFaces, i)
265 label cell0 = layer*surfaceFaces.size() + i;
266 label cell1 = (layer+1)*surfaceFaces.size() + i;
268 eCells[cell0][nCellFaces[cell0]++] = facei;
269 eCells[cell1][nCellFaces[cell1]++] = facei;
276 // External side faces
277 for (label layer=0; layer<nLayers; layer++)
279 // Side faces across layer
280 for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
282 // Get patch faces using edge
283 const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
285 // Get cells on this layer
286 label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
288 eCells[cell0][nCellFaces[cell0]++] = facei;
295 forAll(surfaceFaces, i)
297 eCells[i][nCellFaces[i]++] = facei;
303 forAll(surfaceFaces, i)
305 label cell0 = (nLayers-1)*surfaceFaces.size() + i;
307 eCells[cell0][nCellFaces[cell0]++] = facei;
312 // return points for transferring
313 return xferMove(eCells);
317 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
322 template<class> class FaceList,
325 Foam::extrudedMesh::extrudedMesh
328 const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
329 const extrudeModel& model
335 extrudedPoints(extrudePatch, model),
336 extrudedFaces(extrudePatch, model),
337 extrudedCells(extrudePatch, model)
341 List<polyPatch*> patches(3);
343 label facei = nInternalFaces();
347 *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
349 patches[0] = new wallPolyPatch
360 patches[1] = new polyPatch
369 facei += extrudePatch.size();
371 patches[2] = new polyPatch
384 // ************************************************************************* //