hexBlock.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-2018 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 "hexBlock.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 label hexBlock::vtxLabel(label a, label b, label c) const
36 {
37  return (a + b*(xDim_ + 1) + c*(xDim_ + 1)*(yDim_ + 1));
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 // Construct from components
44 hexBlock::hexBlock(const label nx, const label ny, const label nz)
45 :
46  xDim_(nx),
47  yDim_(ny),
48  zDim_(nz),
49  blockHandedness_(noPoints),
50  points_((xDim_ + 1)*(yDim_ + 1)*(zDim_ + 1))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 void hexBlock::readPoints(Istream& is)
57 {
58  forAll(points_, i)
59  {
60  is >> points_[i].x() >> points_[i].y() >> points_[i].z();
61  }
62 
63  // Calculate the handedness of the block
64  vector i = points_[xDim_] - points_[0];
65  vector j = points_[(xDim_ + 1)*yDim_] - points_[0];
66  vector k = points_[(xDim_ + 1)*(yDim_ + 1)*zDim_] - points_[0];
67 
68  if (((i ^ j) & k) > 0)
69  {
70  Info<< "right-handed block" << endl;
71  blockHandedness_ = right;
72  }
73  else
74  {
75  Info<< "left-handed block" << endl;
76  blockHandedness_ = left;
77  }
78 }
79 
80 
82 {
83  labelListList result(xDim_*yDim_*zDim_);
84 
85  label cellNo = 0;
86 
87  if (blockHandedness_ == right)
88  {
89  for (label k = 0; k <= zDim_ - 1; k++)
90  {
91  for (label j = 0; j <= yDim_ - 1; j++)
92  {
93  for (label i = 0; i <= xDim_ - 1; i++)
94  {
95  labelList& hexLabels = result[cellNo];
96  hexLabels.setSize(8);
97 
98  hexLabels[0] = vtxLabel(i, j, k);
99  hexLabels[1] = vtxLabel(i+1, j, k);
100  hexLabels[2] = vtxLabel(i+1, j+1, k);
101  hexLabels[3] = vtxLabel(i, j+1, k);
102  hexLabels[4] = vtxLabel(i, j, k+1);
103  hexLabels[5] = vtxLabel(i+1, j, k+1);
104  hexLabels[6] = vtxLabel(i+1, j+1, k+1);
105  hexLabels[7] = vtxLabel(i, j+1, k+1);
106 
107  cellNo++;
108  }
109  }
110  }
111  }
112  else if (blockHandedness_ == left)
113  {
114  for (label k = 0; k <= zDim_ - 1; k++)
115  {
116  for (label j = 0; j <= yDim_ - 1; j++)
117  {
118  for (label i = 0; i <= xDim_ - 1; i++)
119  {
120  labelList& hexLabels = result[cellNo];
121  hexLabels.setSize(8);
122 
123  hexLabels[0] = vtxLabel(i, j, k+1);
124  hexLabels[1] = vtxLabel(i+1, j, k+1);
125  hexLabels[2] = vtxLabel(i+1, j+1, k+1);
126  hexLabels[3] = vtxLabel(i, j+1, k+1);
127  hexLabels[4] = vtxLabel(i, j, k);
128  hexLabels[5] = vtxLabel(i+1, j, k);
129  hexLabels[6] = vtxLabel(i+1, j+1, k);
130  hexLabels[7] = vtxLabel(i, j+1, k);
131 
132  cellNo++;
133  }
134  }
135  }
136  }
137  else
138  {
140  << "Unable to determine block handedness as points "
141  << "have not been read in yet"
142  << abort(FatalError);
143  }
144 
145  return result;
146 }
147 
148 
149 // Return block patch faces given direction and range limits
150 // From the cfx manual: direction
151 // 0 = solid (3-D patch),
152 // 1 = high i, 2 = high j, 3 = high k
153 // 4 = low i, 5 = low j, 6 = low k
154 faceList hexBlock::patchFaces(const label direc, const labelList& range) const
155 {
156  if (range.size() != 6)
157  {
159  << "Invalid size of the range array: " << range.size()
160  << ". Should be 6 (xMin, xMax, yMin, yMax, zMin, zMax"
161  << abort(FatalError);
162  }
163 
164  label xMinRange = range[0];
165  label xMaxRange = range[1];
166  label yMinRange = range[2];
167  label yMaxRange = range[3];
168  label zMinRange = range[4];
169  label zMaxRange = range[5];
170 
171  faceList result(0);
172 
173  switch (direc)
174  {
175  case 1:
176  {
177  // high i = xmax
178 
179  result.setSize
180  (
181  (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
182  );
183 
184  label p = 0;
185  for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
186  {
187  for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
188  {
189  result[p].setSize(4);
190 
191  // set the points
192  result[p][0] = vtxLabel(xDim_, j, k);
193  result[p][1] = vtxLabel(xDim_, j+1, k);
194  result[p][2] = vtxLabel(xDim_, j+1, k+1);
195  result[p][3] = vtxLabel(xDim_, j, k+1);
196 
197  p++;
198  }
199  }
200 
201  result.setSize(p);
202  break;
203  }
204 
205  case 2:
206  {
207  // high j = ymax
208  result.setSize
209  (
210  (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
211  );
212 
213  label p = 0;
214  for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
215  {
216  for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
217  {
218  result[p].setSize(4);
219 
220  // set the points
221  result[p][0] = vtxLabel(i, yDim_, k);
222  result[p][1] = vtxLabel(i, yDim_, k + 1);
223  result[p][2] = vtxLabel(i + 1, yDim_, k + 1);
224  result[p][3] = vtxLabel(i + 1, yDim_, k);
225 
226  p++;
227  }
228  }
229 
230  result.setSize(p);
231  break;
232  }
233 
234  case 3:
235  {
236  // high k = zmax
237  result.setSize
238  (
239  (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
240  );
241 
242  label p = 0;
243  for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
244  {
245  for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
246  {
247  result[p].setSize(4);
248 
249  // set the points
250  result[p][0] = vtxLabel(i, j, zDim_);
251  result[p][1] = vtxLabel(i + 1, j, zDim_);
252  result[p][2] = vtxLabel(i + 1, j + 1, zDim_);
253  result[p][3] = vtxLabel(i, j + 1, zDim_);
254 
255  p++;
256  }
257  }
258 
259  result.setSize(p);
260  break;
261  }
262 
263  case 4:
264  {
265  // low i = xmin
266  result.setSize
267  (
268  (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
269  );
270 
271  label p = 0;
272  for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
273  {
274  for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
275  {
276  result[p].setSize(4);
277 
278  // set the points
279  result[p][0] = vtxLabel(0, j, k);
280  result[p][1] = vtxLabel(0, j, k + 1);
281  result[p][2] = vtxLabel(0, j + 1, k + 1);
282  result[p][3] = vtxLabel(0, j + 1, k);
283 
284  p++;
285  }
286  }
287 
288  result.setSize(p);
289  break;
290  }
291 
292  case 5:
293  {
294  // low j = ymin
295  result.setSize
296  (
297  (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
298  );
299 
300  label p = 0;
301  for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
302  {
303  for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
304  {
305  result[p].setSize(4);
306 
307  // set the points
308  result[p][0] = vtxLabel(i, 0, k);
309  result[p][1] = vtxLabel(i + 1, 0, k);
310  result[p][2] = vtxLabel(i + 1, 0, k + 1);
311  result[p][3] = vtxLabel(i, 0, k + 1);
312 
313  p++;
314  }
315  }
316 
317  result.setSize(p);
318  break;
319  }
320 
321  case 6:
322  {
323  // low k = zmin
324  result.setSize
325  (
326  (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
327  );
328 
329  label p = 0;
330  for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
331  {
332  for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
333  {
334  result[p].setSize(4);
335 
336  // set the points
337  result[p][0] = vtxLabel(i, j, 0);
338  result[p][1] = vtxLabel(i, j + 1, 0);
339  result[p][2] = vtxLabel(i + 1, j + 1, 0);
340  result[p][3] = vtxLabel(i + 1, j, 0);
341 
342  p++;
343  }
344  }
345 
346  result.setSize(p);
347  break;
348  }
349 
350  default:
351  {
353  << "direction out of range (1 to 6): " << direc
354  << abort(FatalError);
355  }
356  }
357 
358  // Correct the face orientation based on the handedness of the block.
359  // Do nothing for the right-handed block
360  if (blockHandedness_ == noPoints)
361  {
363  << "Unable to determine block handedness as points "
364  << "have not been read in yet"
365  << abort(FatalError);
366  }
367  else if (blockHandedness_ == left)
368  {
369  // turn all faces inside out
370  forAll(result, facei)
371  {
372  result[facei].flip();
373  }
374  }
375 
376  return result;
377 }
378 
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 } // End namespace Foam
383 
384 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label k
Boltzmann constant.
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
faceList patchFaces(label direc, const labelList &range) const
Return block patch faces given direction and range limits.
void setSize(const label)
Reset size of List.
Definition: List.C:281
messageStream Info
void readPoints(Istream &)
Read block points.
volScalarField & p
hexBlock(const label nx, const label ny, const label nz)
Construct from components.
Namespace for OpenFOAM.
labelListList blockCells() const
Return block cells.