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