1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
27 \*---------------------------------------------------------------------------*/
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 label hexBlock::vtxLabel(label a, label b, label c) const
40 return (a + b*(xDim_ + 1) + c*(xDim_ + 1)*(yDim_ + 1));
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 // Construct from components
47 hexBlock::hexBlock(const label nx, const label ny, const label nz)
52 blockHandedness_(noPoints),
53 points_((xDim_ + 1)*(yDim_ + 1)*(zDim_ + 1))
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 void hexBlock::readPoints(Istream& is)
63 is >> points_[i].x() >> points_[i].y() >> points_[i].z();
66 // Calculate the handedness of the block
67 vector i = points_[xDim_] - points_[0];
68 vector j = points_[(xDim_ + 1)*yDim_] - points_[0];
69 vector k = points_[(xDim_ + 1)*(yDim_ + 1)*zDim_] - points_[0];
71 if (((i ^ j) & k) > 0)
73 Info << "right-handed block" << endl;
74 blockHandedness_ = right;
78 Info << "left-handed block" << endl;
79 blockHandedness_ = left;
84 labelListList hexBlock::blockCells() const
86 labelListList result(xDim_*yDim_*zDim_);
90 if (blockHandedness_ == right)
92 for (label k = 0; k <= zDim_ - 1; k++)
94 for (label j = 0; j <= yDim_ - 1; j++)
96 for (label i = 0; i <= xDim_ - 1; i++)
98 labelList& hexLabels = result[cellNo];
101 hexLabels[0] = vtxLabel(i, j, k);
102 hexLabels[1] = vtxLabel(i+1, j, k);
103 hexLabels[2] = vtxLabel(i+1, j+1, k);
104 hexLabels[3] = vtxLabel(i, j+1, k);
105 hexLabels[4] = vtxLabel(i, j, k+1);
106 hexLabels[5] = vtxLabel(i+1, j, k+1);
107 hexLabels[6] = vtxLabel(i+1, j+1, k+1);
108 hexLabels[7] = vtxLabel(i, j+1, k+1);
115 else if (blockHandedness_ == left)
117 for (label k = 0; k <= zDim_ - 1; k++)
119 for (label j = 0; j <= yDim_ - 1; j++)
121 for (label i = 0; i <= xDim_ - 1; i++)
123 labelList& hexLabels = result[cellNo];
124 hexLabels.setSize(8);
126 hexLabels[0] = vtxLabel(i, j, k+1);
127 hexLabels[1] = vtxLabel(i+1, j, k+1);
128 hexLabels[2] = vtxLabel(i+1, j+1, k+1);
129 hexLabels[3] = vtxLabel(i, j+1, k+1);
130 hexLabels[4] = vtxLabel(i, j, k);
131 hexLabels[5] = vtxLabel(i+1, j, k);
132 hexLabels[6] = vtxLabel(i+1, j+1, k);
133 hexLabels[7] = vtxLabel(i, j+1, k);
142 FatalErrorIn("hexBlock::cellShapes()")
143 << "Unable to determine block handedness as points "
144 << "have not been read in yet"
145 << abort(FatalError);
152 // Return block patch faces given direction and range limits
153 // From the cfx manual: direction
154 // 0 = solid (3-D patch),
155 // 1 = high i, 2 = high j, 3 = high k
156 // 4 = low i, 5 = low j, 6 = low k
157 faceList hexBlock::patchFaces(const label direc, const labelList& range) const
159 if (range.size() != 6)
163 "patchFaces(const label direc, const labelList& range) const"
164 ) << "Invalid size of the range array: " << range.size()
165 << ". Should be 6 (xMin, xMax, yMin, yMax, zMin, zMax"
166 << abort(FatalError);
169 label xMinRange = range[0];
170 label xMaxRange = range[1];
171 label yMinRange = range[2];
172 label yMaxRange = range[3];
173 label zMinRange = range[4];
174 label zMaxRange = range[5];
186 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
190 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
192 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
194 result[p].setSize(4);
197 result[p][0] = vtxLabel(xDim_, j, k);
198 result[p][1] = vtxLabel(xDim_, j+1, k);
199 result[p][2] = vtxLabel(xDim_, j+1, k+1);
200 result[p][3] = vtxLabel(xDim_, j, k+1);
215 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
219 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
221 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
223 result[p].setSize(4);
226 result[p][0] = vtxLabel(i, yDim_, k);
227 result[p][1] = vtxLabel(i, yDim_, k + 1);
228 result[p][2] = vtxLabel(i + 1, yDim_, k + 1);
229 result[p][3] = vtxLabel(i + 1, yDim_, k);
244 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
248 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
250 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
252 result[p].setSize(4);
255 result[p][0] = vtxLabel(i, j, zDim_);
256 result[p][1] = vtxLabel(i + 1, j, zDim_);
257 result[p][2] = vtxLabel(i + 1, j + 1, zDim_);
258 result[p][3] = vtxLabel(i, j + 1, zDim_);
273 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
277 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
279 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
281 result[p].setSize(4);
284 result[p][0] = vtxLabel(0, j, k);
285 result[p][1] = vtxLabel(0, j, k + 1);
286 result[p][2] = vtxLabel(0, j + 1, k + 1);
287 result[p][3] = vtxLabel(0, j + 1, k);
302 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
306 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
308 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
310 result[p].setSize(4);
313 result[p][0] = vtxLabel(i, 0, k);
314 result[p][1] = vtxLabel(i + 1, 0, k);
315 result[p][2] = vtxLabel(i + 1, 0, k + 1);
316 result[p][3] = vtxLabel(i, 0, k + 1);
331 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
335 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
337 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
339 result[p].setSize(4);
342 result[p][0] = vtxLabel(i, j, 0);
343 result[p][1] = vtxLabel(i, j + 1, 0);
344 result[p][2] = vtxLabel(i + 1, j + 1, 0);
345 result[p][3] = vtxLabel(i + 1, j, 0);
359 "patchFaces(const label direc, const labelList& range) const"
360 ) << "direction out of range (1 to 6): " << direc
361 << abort(FatalError);
365 // Correct the face orientation based on the handedness of the block.
366 // Do nothing for the right-handed block
367 if (blockHandedness_ == noPoints)
371 "patchFaces(const label direc, const labelList& range) const"
372 ) << "Unable to determine block handedness as points "
373 << "have not been read in yet"
374 << abort(FatalError);
376 else if (blockHandedness_ == left)
378 // turn all faces inside out
379 forAll (result, faceI)
381 result[faceI] = result[faceI].reverseFace();
389 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 } // End namespace Foam
393 // ************************************************************************* //