ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / mesh / conversion / cfx4ToFoam / hexBlock.C
blob1bbf6f5569224bbac30e59b5c94a744c0d6bbe3c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 \*---------------------------------------------------------------------------*/
26 #include "hexBlock.H"
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 namespace Foam
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)
46     xDim_(nx),
47     yDim_(ny),
48     zDim_(nz),
49     blockHandedness_(noPoints),
50     points_((xDim_ + 1)*(yDim_ + 1)*(zDim_ + 1))
54 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
56 void hexBlock::readPoints(Istream& is)
58     forAll(points_, i)
59     {
60         is  >> points_[i].x() >> points_[i].y() >> points_[i].z();
61     }
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)
69     {
70         Info<< "right-handed block" << endl;
71         blockHandedness_ = right;
72     }
73     else
74     {
75         Info<< "left-handed block" << endl;
76         blockHandedness_ = left;
77     }
81 labelListList hexBlock::blockCells() const
83     labelListList result(xDim_*yDim_*zDim_);
85     label cellNo = 0;
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);
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);
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);
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);
132                     cellNo++;
133                 }
134             }
135         }
136     }
137     else
138     {
139         FatalErrorIn("hexBlock::cellShapes()")
140             << "Unable to determine block handedness as points "
141             << "have not been read in yet"
142             << abort(FatalError);
143     }
145     return result;
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)
157     {
158         FatalErrorIn
159         (
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);
164     }
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];
173     faceList result(0);
175     switch (direc)
176     {
177         case 1:
178         {
179             // high i = xmax
181             result.setSize
182             (
183                 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
184             );
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);
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);
199                     p++;
200                 }
201             }
203             result.setSize(p);
204             break;
205         }
207         case 2:
208         {
209             // high j = ymax
210             result.setSize
211             (
212                 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
213             );
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);
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);
228                     p++;
229                 }
230             }
232             result.setSize(p);
233             break;
234         }
236         case 3:
237         {
238             // high k = zmax
239             result.setSize
240             (
241                 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
242             );
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);
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_);
257                     p++;
258                 }
259             }
261             result.setSize(p);
262             break;
263         }
265         case 4:
266         {
267             // low i = xmin
268             result.setSize
269             (
270                 (yMaxRange - yMinRange + 1)*(zMaxRange - zMinRange + 1)
271             );
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);
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);
286                     p++;
287                 }
288             }
290             result.setSize(p);
291             break;
292         }
294         case 5:
295         {
296             // low j = ymin
297             result.setSize
298             (
299                 (xMaxRange - xMinRange + 1)*(zMaxRange - zMinRange + 1)
300             );
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);
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);
315                     p++;
316                 }
317             }
319             result.setSize(p);
320             break;
321         }
323         case 6:
324         {
325             // low k = zmin
326             result.setSize
327             (
328                 (xMaxRange - xMinRange + 1)*(yMaxRange - yMinRange + 1)
329             );
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);
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);
344                     p++;
345                 }
346             }
348             result.setSize(p);
349             break;
350         }
352         default:
353         {
354             FatalErrorIn
355             (
356                 "patchFaces(const label direc, const labelList& range) const"
357             )   << "direction out of range (1 to 6): " << direc
358                 << abort(FatalError);
359         }
360     }
362     // Correct the face orientation based on the handedness of the block.
363     // Do nothing for the right-handed block
364     if (blockHandedness_ == noPoints)
365     {
366         FatalErrorIn
367         (
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);
372     }
373     else if (blockHandedness_ == left)
374     {
375         // turn all faces inside out
376         forAll(result, faceI)
377         {
378             result[faceI].flip();
379         }
380     }
382     return result;
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 } // End namespace Foam
390 // ************************************************************************* //