extrudedMeshTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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<class FaceList, class PointField>
32 Foam::pointField Foam::extrudedMesh::extrudedPoints
33 (
34  const PrimitivePatch<FaceList, PointField>& extrudePatch,
35  const extrudeModel& model
36 )
37 {
38  const pointField& surfacePoints = extrudePatch.localPoints();
39  const vectorField& surfaceNormals = extrudePatch.pointNormals();
40 
41  const label nLayers = model.nLayers();
42 
43  pointField ePoints((nLayers + 1)*surfacePoints.size());
44 
45  for (label layer=0; layer<=nLayers; layer++)
46  {
47  label offset = layer*surfacePoints.size();
48 
49  forAll(surfacePoints, i)
50  {
51  ePoints[offset + i] = model
52  (
53  surfacePoints[i],
54  surfaceNormals[i],
55  layer
56  );
57  }
58  }
59 
60  return ePoints;
61 }
62 
63 
64 template<class FaceList, class PointField>
65 Foam::faceList Foam::extrudedMesh::extrudedFaces
66 (
67  const PrimitivePatch<FaceList, PointField>& extrudePatch,
68  const extrudeModel& model
69 )
70 {
71  const pointField& surfacePoints = extrudePatch.localPoints();
72  const List<face>& surfaceFaces = extrudePatch.localFaces();
73  const edgeList& surfaceEdges = extrudePatch.edges();
74  const label nInternalEdges = extrudePatch.nInternalEdges();
75 
76  const label nLayers = model.nLayers();
77 
78  label nFaces =
79  (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
80 
81  faceList eFaces(nFaces);
82 
83  labelList quad(4);
84  label facei = 0;
85 
86  // Internal faces
87  for (label layer=0; layer<nLayers; layer++)
88  {
89  label currentLayerOffset = layer*surfacePoints.size();
90  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
91 
92  // Vertical faces from layer to layer+1
93  for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
94  {
95  const edge& e = surfaceEdges[edgeI];
96  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
97 
98  face& f = eFaces[facei++];
99  f.setSize(4);
100 
101  if
102  (
103  (edgeFaces[0] < edgeFaces[1])
104  == sameOrder(surfaceFaces[edgeFaces[0]], e)
105  )
106  {
107  f[0] = e[0] + currentLayerOffset;
108  f[1] = e[1] + currentLayerOffset;
109  f[2] = e[1] + nextLayerOffset;
110  f[3] = e[0] + nextLayerOffset;
111  }
112  else
113  {
114  f[0] = e[1] + currentLayerOffset;
115  f[1] = e[0] + currentLayerOffset;
116  f[2] = e[0] + nextLayerOffset;
117  f[3] = e[1] + nextLayerOffset;
118  }
119  }
120 
121  // Faces between layer and layer+1
122  if (layer < nLayers-1)
123  {
124  forAll(surfaceFaces, i)
125  {
126  eFaces[facei++] =
127  face
128  (
129  surfaceFaces[i] //.reverseFace()
130  + nextLayerOffset
131  );
132  }
133  }
134  }
135 
136  // External side faces
137  for (label layer=0; layer<nLayers; layer++)
138  {
139  label currentLayerOffset = layer*surfacePoints.size();
140  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
141 
142  // Side faces across layer
143  for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
144  {
145  const edge& e = surfaceEdges[edgeI];
146  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
147 
148  face& f = eFaces[facei++];
149  f.setSize(4);
150 
151  if (sameOrder(surfaceFaces[edgeFaces[0]], e))
152  {
153  f[0] = e[0] + currentLayerOffset;
154  f[1] = e[1] + currentLayerOffset;
155  f[2] = e[1] + nextLayerOffset;
156  f[3] = e[0] + nextLayerOffset;
157  }
158  else
159  {
160  f[0] = e[1] + currentLayerOffset;
161  f[1] = e[0] + currentLayerOffset;
162  f[2] = e[0] + nextLayerOffset;
163  f[3] = e[1] + nextLayerOffset;
164  }
165  }
166  }
167 
168  // Bottom faces
169  forAll(surfaceFaces, i)
170  {
171  eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
172  }
173 
174  // Top faces
175  forAll(surfaceFaces, i)
176  {
177  eFaces[facei++] =
178  face
179  (
180  surfaceFaces[i]
181  + nLayers*surfacePoints.size()
182  );
183  }
184 
185  return eFaces;
186 }
187 
188 
189 template<class FaceList, class PointField>
190 Foam::cellList Foam::extrudedMesh::extrudedCells
191 (
192  const PrimitivePatch<FaceList, PointField>& extrudePatch,
193  const extrudeModel& model
194 )
195 {
196  const List<face>& surfaceFaces = extrudePatch.localFaces();
197  const edgeList& surfaceEdges = extrudePatch.edges();
198  const label nInternalEdges = extrudePatch.nInternalEdges();
199 
200  const label nLayers = model.nLayers();
201 
202  cellList eCells(nLayers*surfaceFaces.size());
203 
204  // Size the cells
205  forAll(surfaceFaces, i)
206  {
207  const face& f = surfaceFaces[i];
208 
209  for (label layer=0; layer<nLayers; layer++)
210  {
211  eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
212  }
213  }
214 
215  // Current face count per cell.
216  labelList nCellFaces(eCells.size(), 0);
217 
218 
219  label facei = 0;
220 
221  for (label layer=0; layer<nLayers; layer++)
222  {
223  // Side faces from layer to layer+1
224  for (label i=0; i<nInternalEdges; i++)
225  {
226  // Get patch faces using edge
227  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
228 
229  // Get cells on this layer
230  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
231  label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
232 
233  eCells[cell0][nCellFaces[cell0]++] = facei;
234  eCells[cell1][nCellFaces[cell1]++] = facei;
235 
236  facei++;
237  }
238 
239  // Faces between layer and layer+1
240  if (layer < nLayers-1)
241  {
242  forAll(surfaceFaces, i)
243  {
244  label cell0 = layer*surfaceFaces.size() + i;
245  label cell1 = (layer+1)*surfaceFaces.size() + i;
246 
247  eCells[cell0][nCellFaces[cell0]++] = facei;
248  eCells[cell1][nCellFaces[cell1]++] = facei;
249 
250  facei++;
251  }
252  }
253  }
254 
255  // External side faces
256  for (label layer=0; layer<nLayers; layer++)
257  {
258  // Side faces across layer
259  for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
260  {
261  // Get patch faces using edge
262  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
263 
264  // Get cells on this layer
265  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
266 
267  eCells[cell0][nCellFaces[cell0]++] = facei;
268 
269  facei++;
270  }
271  }
272 
273  // Top faces
274  forAll(surfaceFaces, i)
275  {
276  eCells[i][nCellFaces[i]++] = facei;
277 
278  facei++;
279  }
280 
281  // Bottom faces
282  forAll(surfaceFaces, i)
283  {
284  label cell0 = (nLayers-1)*surfaceFaces.size() + i;
285 
286  eCells[cell0][nCellFaces[cell0]++] = facei;
287 
288  facei++;
289  }
290 
291  return eCells;
292 }
293 
294 
295 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
296 
297 template<class FaceList, class PointField>
299 (
300  const IOobject& io,
301  const PrimitivePatch<FaceList, PointField>& extrudePatch,
302  const extrudeModel& model
303 )
304 :
305  polyMesh
306  (
307  io,
308  extrudedPoints(extrudePatch, model),
309  extrudedFaces(extrudePatch, model),
310  extrudedCells(extrudePatch, model)
311  ),
312  model_(model)
313 {
314  List<polyPatch*> patches(3);
315 
316  label facei = nInternalFaces();
317 
318  label sz =
319  model_.nLayers()
320  *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
321 
322  patches[0] = new wallPolyPatch
323  (
324  "sides",
325  sz,
326  facei,
327  0,
328  boundaryMesh(),
329  wallPolyPatch::typeName
330  );
331 
332  facei += sz;
333 
334  patches[1] = new polyPatch
335  (
336  "originalPatch",
337  extrudePatch.size(),
338  facei,
339  1,
340  boundaryMesh(),
341  polyPatch::typeName
342  );
343 
344  facei += extrudePatch.size();
345 
346  patches[2] = new polyPatch
347  (
348  "otherSide",
349  extrudePatch.size(),
350  facei,
351  2,
352  boundaryMesh(),
353  polyPatch::typeName
354  );
355 
357 }
358 
359 
360 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
label nInternalFaces() const
extrudedMesh(const IOobject &, const PrimitivePatch< FaceList, PointField > &extrudePatch, const extrudeModel &)
Construct from the primitivePatch to extrude.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
List< face > faceList
Definition: faceListFwd.H:43
patches[0]
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< edge > edgeList
Definition: edgeList.H:38
label nLayers() const
Definition: extrudeModel.C:59
List< label > labelList
A List of labels.
Definition: labelList.H:56
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:909
void setSize(const label)
Reset size of List.
Definition: List.C:281
Field< vector > vectorField
Specialisation of Field<T> for vector.
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
List< cell > cellList
list of cells
Definition: cellList.H:42
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163