extrudedMeshTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "extrudedMesh.H"
27 #include "wallPolyPatch.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template
32 <
33  class Face,
34  template<class> class FaceList,
35  class PointField
36 >
38 (
40  const extrudeModel& model
41 )
42 {
43  const pointField& surfacePoints = extrudePatch.localPoints();
44  const vectorField& surfaceNormals = extrudePatch.pointNormals();
45 
46  const label nLayers = model.nLayers();
47 
48  pointField ePoints((nLayers + 1)*surfacePoints.size());
49 
50  for (label layer=0; layer<=nLayers; layer++)
51  {
52  label offset = layer*surfacePoints.size();
53 
54  forAll(surfacePoints, i)
55  {
56  ePoints[offset + i] = model
57  (
58  surfacePoints[i],
59  surfaceNormals[i],
60  layer
61  );
62  }
63  }
64 
65  // return points for transferring
66  return xferMove(ePoints);
67 }
68 
69 
70 template<class Face, template<class> class FaceList, class PointField>
72 (
74  const extrudeModel& model
75 )
76 {
77  const pointField& surfacePoints = extrudePatch.localPoints();
78  const List<face>& surfaceFaces = extrudePatch.localFaces();
79  const edgeList& surfaceEdges = extrudePatch.edges();
80  const label nInternalEdges = extrudePatch.nInternalEdges();
81 
82  const label nLayers = model.nLayers();
83 
84  label nFaces =
85  (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
86 
87  faceList eFaces(nFaces);
88 
89  labelList quad(4);
90  label facei = 0;
91 
92  // Internal faces
93  for (label layer=0; layer<nLayers; layer++)
94  {
95  label currentLayerOffset = layer*surfacePoints.size();
96  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
97 
98  // Vertical faces from layer to layer+1
99  for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
100  {
101  const edge& e = surfaceEdges[edgeI];
102  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
103 
104  face& f = eFaces[facei++];
105  f.setSize(4);
106 
107  if
108  (
109  (edgeFaces[0] < edgeFaces[1])
110  == sameOrder(surfaceFaces[edgeFaces[0]], e)
111  )
112  {
113  f[0] = e[0] + currentLayerOffset;
114  f[1] = e[1] + currentLayerOffset;
115  f[2] = e[1] + nextLayerOffset;
116  f[3] = e[0] + nextLayerOffset;
117  }
118  else
119  {
120  f[0] = e[1] + currentLayerOffset;
121  f[1] = e[0] + currentLayerOffset;
122  f[2] = e[0] + nextLayerOffset;
123  f[3] = e[1] + nextLayerOffset;
124  }
125  }
126 
127  // Faces between layer and layer+1
128  if (layer < nLayers-1)
129  {
130  forAll(surfaceFaces, i)
131  {
132  eFaces[facei++] =
133  face
134  (
135  surfaceFaces[i] //.reverseFace()
136  + nextLayerOffset
137  );
138  }
139  }
140  }
141 
142  // External side faces
143  for (label layer=0; layer<nLayers; layer++)
144  {
145  label currentLayerOffset = layer*surfacePoints.size();
146  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
147 
148  // Side faces across layer
149  for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
150  {
151  const edge& e = surfaceEdges[edgeI];
152  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
153 
154  face& f = eFaces[facei++];
155  f.setSize(4);
156 
157  if (sameOrder(surfaceFaces[edgeFaces[0]], e))
158  {
159  f[0] = e[0] + currentLayerOffset;
160  f[1] = e[1] + currentLayerOffset;
161  f[2] = e[1] + nextLayerOffset;
162  f[3] = e[0] + nextLayerOffset;
163  }
164  else
165  {
166  f[0] = e[1] + currentLayerOffset;
167  f[1] = e[0] + currentLayerOffset;
168  f[2] = e[0] + nextLayerOffset;
169  f[3] = e[1] + nextLayerOffset;
170  }
171  }
172  }
173 
174  // Bottom faces
175  forAll(surfaceFaces, i)
176  {
177  eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
178  }
179 
180  // Top faces
181  forAll(surfaceFaces, i)
182  {
183  eFaces[facei++] =
184  face
185  (
186  surfaceFaces[i]
187  + nLayers*surfacePoints.size()
188  );
189  }
190 
191  // return points for transferring
192  return xferMove(eFaces);
193 }
194 
195 
196 template<class Face, template<class> class FaceList, class PointField>
198 (
199  const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
200  const extrudeModel& model
201 )
202 {
203  const List<face>& surfaceFaces = extrudePatch.localFaces();
204  const edgeList& surfaceEdges = extrudePatch.edges();
205  const label nInternalEdges = extrudePatch.nInternalEdges();
206 
207  const label nLayers = model.nLayers();
208 
209  cellList eCells(nLayers*surfaceFaces.size());
210 
211  // Size the cells
212  forAll(surfaceFaces, i)
213  {
214  const face& f = surfaceFaces[i];
215 
216  for (label layer=0; layer<nLayers; layer++)
217  {
218  eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
219  }
220  }
221 
222  // Current face count per cell.
223  labelList nCellFaces(eCells.size(), 0);
224 
225 
226  label facei = 0;
227 
228  for (label layer=0; layer<nLayers; layer++)
229  {
230  // Side faces from layer to layer+1
231  for (label i=0; i<nInternalEdges; i++)
232  {
233  // Get patch faces using edge
234  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
235 
236  // Get cells on this layer
237  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
238  label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
239 
240  eCells[cell0][nCellFaces[cell0]++] = facei;
241  eCells[cell1][nCellFaces[cell1]++] = facei;
242 
243  facei++;
244  }
245 
246  // Faces between layer and layer+1
247  if (layer < nLayers-1)
248  {
249  forAll(surfaceFaces, i)
250  {
251  label cell0 = layer*surfaceFaces.size() + i;
252  label cell1 = (layer+1)*surfaceFaces.size() + i;
253 
254  eCells[cell0][nCellFaces[cell0]++] = facei;
255  eCells[cell1][nCellFaces[cell1]++] = facei;
256 
257  facei++;
258  }
259  }
260  }
261 
262  // External side faces
263  for (label layer=0; layer<nLayers; layer++)
264  {
265  // Side faces across layer
266  for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
267  {
268  // Get patch faces using edge
269  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
270 
271  // Get cells on this layer
272  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
273 
274  eCells[cell0][nCellFaces[cell0]++] = facei;
275 
276  facei++;
277  }
278  }
279 
280  // Top faces
281  forAll(surfaceFaces, i)
282  {
283  eCells[i][nCellFaces[i]++] = facei;
284 
285  facei++;
286  }
287 
288  // Bottom faces
289  forAll(surfaceFaces, i)
290  {
291  label cell0 = (nLayers-1)*surfaceFaces.size() + i;
292 
293  eCells[cell0][nCellFaces[cell0]++] = facei;
294 
295  facei++;
296  }
297 
298  // return points for transferring
299  return xferMove(eCells);
300 }
301 
302 
303 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
304 
305 template
306 <
307  class Face,
308  template<class> class FaceList,
309  class PointField
310 >
312 (
313  const IOobject& io,
314  const PrimitivePatch<Face, FaceList, PointField>& extrudePatch,
315  const extrudeModel& model
316 )
317 :
318  polyMesh
319  (
320  io,
321  extrudedPoints(extrudePatch, model),
322  extrudedFaces(extrudePatch, model),
323  extrudedCells(extrudePatch, model)
324  ),
325  model_(model)
326 {
328 
329  label facei = nInternalFaces();
330 
331  label sz =
332  model_.nLayers()
333  *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
334 
335  patches[0] = new wallPolyPatch
336  (
337  "sides",
338  sz,
339  facei,
340  0,
341  boundaryMesh(),
342  wallPolyPatch::typeName
343  );
344 
345  facei += sz;
346 
347  patches[1] = new polyPatch
348  (
349  "originalPatch",
350  extrudePatch.size(),
351  facei,
352  1,
353  boundaryMesh(),
354  polyPatch::typeName
355  );
356 
357  facei += extrudePatch.size();
358 
359  patches[2] = new polyPatch
360  (
361  "otherSide",
362  extrudePatch.size(),
363  facei,
364  2,
365  boundaryMesh(),
366  polyPatch::typeName
367  );
368 
369  addPatches(patches);
370 }
371 
372 
373 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::face::reverseFace
face reverseFace() const
Return face with reverse direction.
Definition: face.C:611
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatchTemplate.C:292
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatchTemplate.C:212
Foam::extrudeModel
Top level extrusion model class.
Definition: extrudeModel.H:51
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatchTemplate.H:299
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
wallPolyPatch.H
Foam::extrudedMesh::extrudedFaces
Xfer< faceList > extrudedFaces(const PrimitivePatch< Face, FaceList, PointField > &extrudePatch, const extrudeModel &)
Construct and return the extruded mesh faces.
Foam::extrudeModel::nLayers
label nLayers() const
Definition: extrudeModel.C:59
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::extrudedMesh::extrudedMesh
extrudedMesh(const extrudedMesh &)
Disallow default bitwise copy construct.
Foam::wallPolyPatch
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Foam::extrudedMesh::extrudedCells
Xfer< cellList > extrudedCells(const PrimitivePatch< Face, FaceList, PointField > &extrudePatch, const extrudeModel &)
Construct and return the extruded mesh cells.
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatchTemplate.C:232
Foam::PrimitivePatch::pointNormals
const Field< PointType > & pointNormals() const
Return point normals for patch.
Definition: PrimitivePatchTemplate.C:540
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
f
labelList f(nPoints)
Foam::List< face >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:59
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
extrudedMesh.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::extrudedMesh::extrudedPoints
Xfer< pointField > extrudedPoints(const PrimitivePatch< Face, FaceList, PointField > &extrudePatch, const extrudeModel &)
Construct and return the extruded mesh points.
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88