1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 label hexBlock::vtxLabel(label a, label b, label c) const
37 return (a + b*(xDim_ + 1) + c*(xDim_ + 1)*(yDim_ + 1));
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 // Construct from components
44 hexBlock::hexBlock(const label nx, const label ny, const label nz)
49 blockHandedness_(noPoints),
50 points_((xDim_ + 1)*(yDim_ + 1)*(zDim_ + 1))
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 void hexBlock::readPoints(Istream& is)
60 is >> points_[i].x() >> points_[i].y() >> points_[i].z();
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];
68 if (((i ^ j) & k) > 0)
70 Info<< "right-handed block" << endl;
71 blockHandedness_ = right;
75 Info<< "left-handed block" << endl;
76 blockHandedness_ = left;
81 labelListList hexBlock::blockCells() const
83 labelListList result(xDim_*yDim_*zDim_);
87 if (blockHandedness_ == right)
89 for (label k = 0; k <= zDim_ - 1; k++)
91 for (label j = 0; j <= yDim_ - 1; j++)
93 for (label i = 0; i <= xDim_ - 1; i++)
95 labelList& hexLabels = result[cellNo];
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);
112 else if (blockHandedness_ == left)
114 for (label k = 0; k <= zDim_ - 1; k++)
116 for (label j = 0; j <= yDim_ - 1; j++)
118 for (label i = 0; i <= xDim_ - 1; i++)
120 labelList& hexLabels = result[cellNo];
121 hexLabels.setSize(8);
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);
139 FatalErrorIn("hexBlock::cellShapes()")
140 << "Unable to determine block handedness as points "
141 << "have not been read in yet"
142 << abort(FatalError);
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
156 if (range.size() != 6)
160 "patchFaces(const label direc, const labelList& range) const"
161 ) << "Invalid size of the range array: " << range.size()
162 << ". Should be 6 (xMin, xMax, yMin, yMax, zMin, zMax"
163 << abort(FatalError);
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];
183 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
187 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
189 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
191 result[p].setSize(4);
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);
212 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
216 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
218 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
220 result[p].setSize(4);
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);
241 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
245 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
247 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
249 result[p].setSize(4);
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_);
270 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
274 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
276 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
278 result[p].setSize(4);
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);
299 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
303 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
305 for (label k = zMinRange - 1; k <= zMaxRange - 1; k++)
307 result[p].setSize(4);
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);
328 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
332 for (label i = xMinRange - 1; i <= xMaxRange - 1; i++)
334 for (label j = yMinRange - 1; j <= yMaxRange - 1; j++)
336 result[p].setSize(4);
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);
356 "patchFaces(const label direc, const labelList& range) const"
357 ) << "direction out of range (1 to 6): " << direc
358 << abort(FatalError);
362 // Correct the face orientation based on the handedness of the block.
363 // Do nothing for the right-handed block
364 if (blockHandedness_ == noPoints)
368 "patchFaces(const label direc, const labelList& range) const"
369 ) << "Unable to determine block handedness as points "
370 << "have not been read in yet"
371 << abort(FatalError);
373 else if (blockHandedness_ == left)
375 // turn all faces inside out
376 forAll(result, faceI)
378 result[faceI].flip();
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 } // End namespace Foam
390 // ************************************************************************* //