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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2012 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "extrudedMesh.H"
30 #include "wallPolyPatch.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class FaceList, class PointField>
35 Foam::pointField Foam::extrudedMesh::extrudedPoints
36 (
37  const PrimitivePatch<FaceList, PointField>& extrudePatch,
38  const extrudeModel& model
39 )
40 {
41  const pointField& surfacePoints = extrudePatch.localPoints();
42  const vectorField& surfaceNormals = extrudePatch.pointNormals();
43 
44  const label nLayers = model.nLayers();
45 
46  pointField ePoints((nLayers + 1)*surfacePoints.size());
47 
48  for (label layer=0; layer <= nLayers; ++layer)
49  {
50  const label offset = layer*surfacePoints.size();
51 
52  forAll(surfacePoints, i)
53  {
54  ePoints[offset + i] = model
55  (
56  surfacePoints[i],
57  surfaceNormals[i],
58  layer
59  );
60  }
61  }
62 
63  return ePoints;
64 }
65 
66 
67 template<class FaceList, class PointField>
68 Foam::faceList Foam::extrudedMesh::extrudedFaces
69 (
70  const PrimitivePatch<FaceList, PointField>& extrudePatch,
71  const extrudeModel& model
72 )
73 {
74  const pointField& surfacePoints = extrudePatch.localPoints();
75  const List<face>& surfaceFaces = extrudePatch.localFaces();
76  const edgeList& surfaceEdges = extrudePatch.edges();
77  const label nInternalEdges = extrudePatch.nInternalEdges();
78 
79  const label nLayers = model.nLayers();
80 
81  label nFaces =
82  (nLayers + 1)*surfaceFaces.size() + nLayers*surfaceEdges.size();
83 
84  faceList eFaces(nFaces);
85 
86  labelList quad(4);
87  label facei = 0;
88 
89  // Internal faces
90  for (label layer=0; layer<nLayers; layer++)
91  {
92  label currentLayerOffset = layer*surfacePoints.size();
93  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
94 
95  // Vertical faces from layer to layer+1
96  for (label edgeI=0; edgeI<nInternalEdges; edgeI++)
97  {
98  const edge& e = surfaceEdges[edgeI];
99  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
100 
101  face& f = eFaces[facei++];
102  f.setSize(4);
103 
104  if
105  (
106  (edgeFaces[0] < edgeFaces[1])
107  == sameOrder(surfaceFaces[edgeFaces[0]], e)
108  )
109  {
110  f[0] = e[0] + currentLayerOffset;
111  f[1] = e[1] + currentLayerOffset;
112  f[2] = e[1] + nextLayerOffset;
113  f[3] = e[0] + nextLayerOffset;
114  }
115  else
116  {
117  f[0] = e[1] + currentLayerOffset;
118  f[1] = e[0] + currentLayerOffset;
119  f[2] = e[0] + nextLayerOffset;
120  f[3] = e[1] + nextLayerOffset;
121  }
122  }
123 
124  // Faces between layer and layer+1
125  if (layer < nLayers-1)
126  {
127  forAll(surfaceFaces, i)
128  {
129  eFaces[facei++] =
130  face
131  (
132  surfaceFaces[i] //.reverseFace()
133  + nextLayerOffset
134  );
135  }
136  }
137  }
138 
139  // External side faces
140  for (label layer=0; layer<nLayers; layer++)
141  {
142  label currentLayerOffset = layer*surfacePoints.size();
143  label nextLayerOffset = currentLayerOffset + surfacePoints.size();
144 
145  // Side faces across layer
146  for (label edgeI=nInternalEdges; edgeI<surfaceEdges.size(); edgeI++)
147  {
148  const edge& e = surfaceEdges[edgeI];
149  const labelList& edgeFaces = extrudePatch.edgeFaces()[edgeI];
150 
151  face& f = eFaces[facei++];
152  f.setSize(4);
153 
154  if (sameOrder(surfaceFaces[edgeFaces[0]], e))
155  {
156  f[0] = e[0] + currentLayerOffset;
157  f[1] = e[1] + currentLayerOffset;
158  f[2] = e[1] + nextLayerOffset;
159  f[3] = e[0] + nextLayerOffset;
160  }
161  else
162  {
163  f[0] = e[1] + currentLayerOffset;
164  f[1] = e[0] + currentLayerOffset;
165  f[2] = e[0] + nextLayerOffset;
166  f[3] = e[1] + nextLayerOffset;
167  }
168  }
169  }
170 
171  // Bottom faces
172  forAll(surfaceFaces, i)
173  {
174  eFaces[facei++] = face(surfaceFaces[i]).reverseFace();
175  }
176 
177  // Top faces
178  forAll(surfaceFaces, i)
179  {
180  eFaces[facei++] =
181  face
182  (
183  surfaceFaces[i]
184  + nLayers*surfacePoints.size()
185  );
186  }
187 
188 
189  return eFaces;
190 }
191 
192 
193 template<class FaceList, class PointField>
194 Foam::cellList Foam::extrudedMesh::extrudedCells
195 (
196  const PrimitivePatch<FaceList, PointField>& extrudePatch,
197  const extrudeModel& model
198 )
199 {
200  const List<face>& surfaceFaces = extrudePatch.localFaces();
201  const edgeList& surfaceEdges = extrudePatch.edges();
202  const label nInternalEdges = extrudePatch.nInternalEdges();
203 
204  const label nLayers = model.nLayers();
205 
206  cellList eCells(nLayers*surfaceFaces.size());
207 
208  // Size the cells
209  forAll(surfaceFaces, i)
210  {
211  const face& f = surfaceFaces[i];
212 
213  for (label layer=0; layer<nLayers; layer++)
214  {
215  eCells[i + layer*surfaceFaces.size()].setSize(f.size() + 2);
216  }
217  }
218 
219  // Current face count per cell.
220  labelList nCellFaces(eCells.size(), Zero);
221 
222 
223  label facei = 0;
224 
225  for (label layer=0; layer<nLayers; layer++)
226  {
227  // Side faces from layer to layer+1
228  for (label i=0; i<nInternalEdges; i++)
229  {
230  // Get patch faces using edge
231  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
232 
233  // Get cells on this layer
234  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
235  label cell1 = layer*surfaceFaces.size() + edgeFaces[1];
236 
237  eCells[cell0][nCellFaces[cell0]++] = facei;
238  eCells[cell1][nCellFaces[cell1]++] = facei;
239 
240  facei++;
241  }
242 
243  // Faces between layer and layer+1
244  if (layer < nLayers-1)
245  {
246  forAll(surfaceFaces, i)
247  {
248  label cell0 = layer*surfaceFaces.size() + i;
249  label cell1 = (layer+1)*surfaceFaces.size() + i;
250 
251  eCells[cell0][nCellFaces[cell0]++] = facei;
252  eCells[cell1][nCellFaces[cell1]++] = facei;
253 
254  facei++;
255  }
256  }
257  }
258 
259  // External side faces
260  for (label layer=0; layer<nLayers; layer++)
261  {
262  // Side faces across layer
263  for (label i=nInternalEdges; i<surfaceEdges.size(); i++)
264  {
265  // Get patch faces using edge
266  const labelList& edgeFaces = extrudePatch.edgeFaces()[i];
267 
268  // Get cells on this layer
269  label cell0 = layer*surfaceFaces.size() + edgeFaces[0];
270 
271  eCells[cell0][nCellFaces[cell0]++] = facei;
272 
273  facei++;
274  }
275  }
276 
277  // Top faces
278  forAll(surfaceFaces, i)
279  {
280  eCells[i][nCellFaces[i]++] = facei;
281 
282  facei++;
283  }
284 
285  // Bottom faces
286  forAll(surfaceFaces, i)
287  {
288  label cell0 = (nLayers-1)*surfaceFaces.size() + i;
289 
290  eCells[cell0][nCellFaces[cell0]++] = facei;
291 
292  facei++;
293  }
294 
295  return eCells;
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
300 
301 template<class FaceList, class PointField>
302 Foam::extrudedMesh::extrudedMesh
303 (
304  const IOobject& io,
305  const PrimitivePatch<FaceList, PointField>& extrudePatch,
306  const extrudeModel& model
307 )
308 :
309  polyMesh
310  (
311  io,
312  extrudedPoints(extrudePatch, model),
313  extrudedFaces(extrudePatch, model),
314  extrudedCells(extrudePatch, model)
315  ),
316  model_(model)
317 {
318  List<polyPatch*> patches(3);
319 
320  label facei = nInternalFaces();
321 
322  label sz =
323  model_.nLayers()
324  *(extrudePatch.nEdges() - extrudePatch.nInternalEdges());
325 
326  patches[0] = new wallPolyPatch
327  (
328  "sides",
329  sz,
330  facei,
331  0,
332  boundaryMesh(),
333  wallPolyPatch::typeName
334  );
335 
336  facei += sz;
337 
338  patches[1] = new polyPatch
339  (
340  "originalPatch",
341  extrudePatch.size(),
342  facei,
343  1,
344  boundaryMesh(),
345  polyPatch::typeName
346  );
347 
348  facei += extrudePatch.size();
349 
350  patches[2] = new polyPatch
351  (
352  "otherSide",
353  extrudePatch.size(),
354  facei,
355  2,
356  boundaryMesh(),
357  polyPatch::typeName
358  );
359 
361 }
362 
363 
364 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
void setSize(const label n)
Alias for resize()
Definition: List.H:289
List< edge > edgeList
A List of edges.
Definition: edgeList.H:60
label nLayers() const
Return the number of layers.
Definition: extrudeModel.C:52
label nInternalFaces() const noexcept
Number of internal faces.
labelList f(nPoints)
const polyBoundaryMesh & patches
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:959
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:41
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157